CASToR  3.2
Tomographic Reconstruction (PET/SPECT/CT)
code/src/image/oImageSpace.cc
Go to the documentation of this file.
1 
8 #include "oImageSpace.hh"
9 #include "vOptimizer.hh"
10 
11 // =====================================================================
12 // ---------------------------------------------------------------------
13 // ---------------------------------------------------------------------
14 // =====================================================================
15 
17 {
18  // Initialize the member variables to their default values
19  m_loadedSensitivity = false;
20  m_loadedMultiModal = false;
21  m_loadedMask = false;
22  mp_ID = NULL;
23  m_verbose = -1;
27  // Set all images pointers to NULL
28  m4p_image = NULL;
29  m4p_forwardImage = NULL;
30  m6p_backwardImage = NULL;
31  m5p_sensitivity = NULL;
32  m4p_attenuation = NULL;
35  m4p_outputImage = NULL;
36  mp_visitedVoxelsImage = NULL;
37  m2p_projectionImage = NULL;
38  m2p_multiModalImage = NULL;
39  mp_maskImage = NULL;
40 }
41 
42 // =====================================================================
43 // ---------------------------------------------------------------------
44 // ---------------------------------------------------------------------
45 // =====================================================================
46 
48 {
49  // Call the deallocation of the miscellaneous images here
51 }
52 
53 // =====================================================================
54 // ---------------------------------------------------------------------
55 // ---------------------------------------------------------------------
56 // =====================================================================
57 
59 {
60  // Verbose
61  if (m_verbose>=3) Cout("oImageSpace::InstantiateImage() -> Initialize to 0."<< endl);
62 
63  // For each temporal dimensions, the number of related images is recovered from the number of
64  // intrinsic time basis functions stored in mp_ID. If no intrinsic temporal basis functions
65  // are initialized (standard reconstruction), this is equivalent to the standard calls to
66  // GetNbTimeFrames() / GetNbRespGates() / GetNbCardGates()
67 
68  // First pointer is for the number of time basis functions
69  m4p_image = (FLTNB****)malloc(mp_ID->GetNbTimeBasisFunctions()*sizeof(FLTNB***));
70  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
71  {
72  // Second pointer is for the number of respiratory basis functions
73  m4p_image[tbf] = (FLTNB***)malloc(mp_ID->GetNbRespBasisFunctions()*sizeof(FLTNB**));
74  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
75  {
76  // Third pointer is for the number of cardiac basis functions
77  m4p_image[tbf][rbf] = (FLTNB**)malloc(mp_ID->GetNbCardBasisFunctions()*sizeof(FLTNB*));
78  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
79  {
80  // Fourth pointer is for the 3D space
81  m4p_image[tbf][rbf][cbf] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
82  // Initialize to 0.
83  for (INTNB v=0; v<mp_ID->GetNbVoxXYZ(); v++) m4p_image[tbf][rbf][cbf][v] = 0.;
84  }
85  }
86  }
87 }
88 
89 // =====================================================================
90 // ---------------------------------------------------------------------
91 // ---------------------------------------------------------------------
92 // =====================================================================
93 
95 {
96  // Check the first pointer
97  if (m4p_image)
98  {
99  // Verbose
100  if (m_verbose>=3) Cout("oImageSpace::DeallocateImage() -> Free memory"<< endl);
101  // Loop on time basis functions
102  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++) if (m4p_image[tbf])
103  {
104  // Loop on respiratory basis functions
105  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++) if (m4p_image[tbf][rbf])
106  {
107  // Loop on cardiac basis functions
108  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++) if (m4p_image[tbf][rbf][cbf])
109  {
110  free(m4p_image[tbf][rbf][cbf]);
111  }
112  free(m4p_image[tbf][rbf]);
113  }
114  free(m4p_image[tbf]);
115  }
116  free(m4p_image);
117  }
118  // Reset the pointer to NULL
119  m4p_image = NULL;
120 }
121 
122 // =====================================================================
123 // ---------------------------------------------------------------------
124 // ---------------------------------------------------------------------
125 // =====================================================================
126 
128 {
129  // Verbose
130  if (m_verbose>=3) Cout("oImageSpace::InstantiateForwardImage() -> Initialize to 0."<< endl);
131 
132  // For each temporal dimensions, the number of related images is recovered from the number of
133  // intrinsic time basis functions stored in mp_ID. If no intrinsic temporal basis functions
134  // are initialized (standard reconstruction), this is equivalent to the standard calls to
135  // GetNbTimeFrames() / GetNbRespGates() / GetNbCardGates()
136 
137  // First pointer is for the number of time basis functions
138  m4p_forwardImage = (FLTNB****)malloc(mp_ID->GetNbTimeBasisFunctions()*sizeof(FLTNB***));
139  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
140  {
141  // Second pointer is for the number of respiratory basis functions
142  m4p_forwardImage[tbf] = (FLTNB***)malloc(mp_ID->GetNbRespBasisFunctions()*sizeof(FLTNB**));
143  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
144  {
145  // Third pointer is for the number of cardiac basis functions
146  m4p_forwardImage[tbf][rbf] = (FLTNB**)malloc(mp_ID->GetNbCardBasisFunctions()*sizeof(FLTNB*));
147  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
148  {
149  // Fourth pointer is for the 3D space
150  m4p_forwardImage[tbf][rbf][cbf] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
151  // Initialize to 0.
152  for (INTNB v=0; v<mp_ID->GetNbVoxXYZ(); v++) m4p_forwardImage[tbf][rbf][cbf][v] = 0.;
153  }
154  }
155  }
156 }
157 
158 // =====================================================================
159 // ---------------------------------------------------------------------
160 // ---------------------------------------------------------------------
161 // =====================================================================
162 
164 {
165  // Check the first pointer
166  if (m4p_forwardImage)
167  {
168  // Verbose
169  if (m_verbose>=3) Cout("oImageSpace::DeallocateForwardImage() -> Free memory"<< endl);
170  // Loop on time basis functions
171  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++) if (m4p_forwardImage[tbf])
172  {
173  // Loop on respiratory basis functions
174  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++) if (m4p_forwardImage[tbf][rbf])
175  {
176  // Loop on cardiac basis functions
177  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++) if (m4p_forwardImage[tbf][rbf][cbf])
178  {
179  free(m4p_forwardImage[tbf][rbf][cbf]);
180  }
181  free(m4p_forwardImage[tbf][rbf]);
182  }
183  free(m4p_forwardImage[tbf]);
184  }
185  free(m4p_forwardImage);
186  }
187 
188  // Reset the pointer to NULL
189  m4p_forwardImage = NULL;
190 }
191 
192 // =====================================================================
193 // ---------------------------------------------------------------------
194 // ---------------------------------------------------------------------
195 // =====================================================================
196 
198 {
199  // Verbose
200  if (m_verbose>=3) Cout("oImageSpace::InstantiateBackwardImageFromDynamicBasis() -> Initialize to 0."<< endl);
201 
202  // For each temporal dimensions, the number of related images is recovered from the number of
203  // intrinsic time basis functions stored in mp_ID. If no intrinsic temporal basis functions
204  // are initialized (standard reconstruction), this is equivalent to the standard calls to
205  // GetNbTimeFrames() / GetNbRespGates() / GetNbCardGates()
206 
207  // First pointer is for the number of backward images (in case the optimizer does need more than 1)
208  m_nbBackwardImages = a_nbBackwardImages;
209  m6p_backwardImage = (FLTNB******)malloc(m_nbBackwardImages*sizeof(FLTNB*****));
210  for (int img=0; img<m_nbBackwardImages; img++)
211  {
212  // Second pointer is for the number of threads, in order to have thread-safe backward projections
213  m6p_backwardImage[img] = (FLTNB*****)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(FLTNB****));
214  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
215  {
216  // The third pointer is for the number of time basis functions
217  m6p_backwardImage[img][th] = (FLTNB****)malloc(mp_ID->GetNbTimeBasisFunctions()*sizeof(FLTNB***));
218  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
219  {
220  // The fourth pointer is for the number of respiratory basis functions
221  m6p_backwardImage[img][th][tbf] = (FLTNB***)malloc(mp_ID->GetNbRespBasisFunctions()*sizeof(FLTNB**));
222  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
223  {
224  // The fifth pointer is for the number of cardiac basis functions
225  m6p_backwardImage[img][th][tbf][rbf] = (FLTNB**)malloc(mp_ID->GetNbCardBasisFunctions()*sizeof(FLTNB*));
226  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
227  {
228  // The sixth pointer is for the 3D space
229  m6p_backwardImage[img][th][tbf][rbf][cbf] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
230  }
231  }
232  }
233  }
234  }
235 }
236 
237 // =====================================================================
238 // ---------------------------------------------------------------------
239 // ---------------------------------------------------------------------
240 // =====================================================================
241 
243 {
244  // Check the first pointer
245  if (m6p_backwardImage)
246  {
247  // Verbose
248  if (m_verbose>=3) Cout("oImageSpace::DeallocateBackwardImageFromDynamicBasis() -> Free memory" << endl);
249  // Loop on backward images (in case the optimizer does need more than 1)
250  for (int img=0; img<m_nbBackwardImages; img++) if (m6p_backwardImage[img]!=NULL)
251  {
252  // Loop on threads
253  for (int th=0 ; th<mp_ID->GetNbThreadsForProjection() ; th++) if (m6p_backwardImage[img][th])
254  {
255  // Loop on time basis functions
256  for (int tbf=0 ; tbf<mp_ID->GetNbTimeBasisFunctions() ; tbf++) if (m6p_backwardImage[img][th][tbf])
257  {
258  // Loop on respiratory basis functions
259  for (int rbf=0 ; rbf<mp_ID->GetNbRespBasisFunctions() ; rbf++) if (m6p_backwardImage[img][th][tbf][rbf])
260  {
261  // Loop on cardiac basis functions
262  for (int cbf=0 ; cbf<mp_ID->GetNbCardBasisFunctions() ; cbf++) if (m6p_backwardImage[img][th][tbf][rbf][cbf])
263  {
264  free(m6p_backwardImage[img][th][tbf][rbf][cbf]);
265  }
266  free(m6p_backwardImage[img][th][tbf][rbf]);
267  }
268  free(m6p_backwardImage[img][th][tbf]);
269  }
270  free(m6p_backwardImage[img][th]);
271  }
272  free(m6p_backwardImage[img]);
273  }
274  free(m6p_backwardImage);
275  }
276  // Reset the pointer to NULL
277  m6p_backwardImage = NULL;
278 }
279 
280 // =====================================================================
281 // ---------------------------------------------------------------------
282 // ---------------------------------------------------------------------
283 // =====================================================================
284 
286 {
287  // Verbose
288  if (m_verbose>=3) Cout("oImageSpace::InstantiateBackwardImageFromDynamicBins() -> Initialize to 0."<< endl);
289  // Not related to an optimizer, so we consider only one image
290  m6p_backwardImage = (FLTNB******)malloc(1*sizeof(FLTNB*****));
291  // Projection threads
292  m6p_backwardImage[0] = (FLTNB*****)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(FLTNB****));
293  for (int th=0 ; th<mp_ID->GetNbThreadsForProjection() ; th++)
294  {
295  // Time frames
296  m6p_backwardImage[0][th] = (FLTNB****)malloc(mp_ID->GetNbTimeFrames()*sizeof(FLTNB***));
297  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
298  {
299  // Respiratory gates
300  m6p_backwardImage[0][th][fr] = (FLTNB***)malloc(mp_ID->GetNb1stMotImgsForLMS(fr)*sizeof(FLTNB**));
301  for (int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
302  {
303  // Cardiac gates
304  m6p_backwardImage[0][th][fr][rg] = (FLTNB**)malloc(mp_ID->GetNb2ndMotImgsForLMS()*sizeof(FLTNB*));
305  for (int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
306  {
307  // Number of voxels
308  m6p_backwardImage[0][th][fr][rg][cg] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
309  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++) m6p_backwardImage[0][th][fr][rg][cg][v] = 0.;
310  }
311  }
312  }
313  }
314 }
315 
316 // =====================================================================
317 // ---------------------------------------------------------------------
318 // ---------------------------------------------------------------------
319 // =====================================================================
320 
322 {
323  // Check the first pointer
324  if (m6p_backwardImage)
325  {
326  // Verbose
327  if (m_verbose>=3) Cout("oImageSpace::DeallocateBackwardImageFromDynamicBins() -> Free memory"<< endl);
328  // We assumed only one image when considering dynamic bins
329  if (m6p_backwardImage[0]!=NULL)
330  {
331  // Loop on threads
332  for (int th=0 ; th<mp_ID->GetNbThreadsForProjection() ; th++) if (m6p_backwardImage[0][th])
333  {
334  // Loop on time frames
335  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++) if (m6p_backwardImage[0][th][fr])
336  {
337  // Loop on respiratory gates
338  for (int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++) if (m6p_backwardImage[0][th][fr][rg])
339  {
340  // Loop on cardiac gates
341  for (int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++) if (m6p_backwardImage[0][th][fr][rg][cg])
342  {
343  free(m6p_backwardImage[0][th][fr][rg][cg]);
344  }
345  free(m6p_backwardImage[0][th][fr][rg]);
346  }
347  free(m6p_backwardImage[0][th][fr]);
348  }
349  free(m6p_backwardImage[0][th]);
350  }
351  free(m6p_backwardImage[0]);
352  }
353  free(m6p_backwardImage);
354  }
355  // Reset the pointer to NULL
356  m6p_backwardImage = NULL;
357 }
358 
359 // =====================================================================
360 // ---------------------------------------------------------------------
361 // ---------------------------------------------------------------------
362 // =====================================================================
363 
364 void oImageSpace::InstantiateSensitivityImage(const string& a_pathToSensitivityImage)
365 {
366  // -----------------------------------------------------------------------------------------------------
367  // Case 1: a path to a sensitivity image is provided, meaning that we will not compute the sensitivity
368  // on-the-fly. In that case, we do not take the multi-threading into account. Also note that in
369  // histogram mode, the same sensitivity image will be used in optimization algorithms for all
370  // subsets.
371  // -----------------------------------------------------------------------------------------------------
372 
373  if (!a_pathToSensitivityImage.empty())
374  {
375  // Verbose
376  if (m_verbose>=3) Cout("oImageSpace::InstantiateSensitivityImage() -> Will be loaded from '" << a_pathToSensitivityImage << "'" << endl);
377  // Allocate for all threads first
378  m5p_sensitivity = (FLTNB*****)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(FLTNB****));
379  // Allocate the rest only for the first thread
380  int thread_0 = 0;
381  m5p_sensitivity[thread_0] = (FLTNB****)malloc(mp_ID->GetNbTimeFrames()*sizeof(FLTNB***));
382  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
383  {
384  m5p_sensitivity[thread_0][fr] = (FLTNB***)malloc(mp_ID->GetNbRespGates()*sizeof(FLTNB**));
385  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
386  {
387  m5p_sensitivity[thread_0][fr][rg] = (FLTNB**)malloc(mp_ID->GetNbCardGates()*sizeof(FLTNB*));
388  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
389  m5p_sensitivity[thread_0][fr][rg][cg] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
390  }
391  }
392  // Make all thread pointers pointing to the first
393  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
394  m5p_sensitivity[th] = m5p_sensitivity[thread_0];
395  // Set the flag for loaded sensitivity to true
396  m_loadedSensitivity = true;
397  }
398 
399  // -----------------------------------------------------------------------------------------------------
400  // Case 2: no sensitivity provided, meaning we will commpute it on-the-fly. We thus need it to be
401  // thread-safe, so we allocate for all threads.
402  // -----------------------------------------------------------------------------------------------------
403 
404  else
405  {
406  // Verbose
407  if (m_verbose>=3) Cout("oImageSpace::InstantiateSensitivityImage() -> For all threads"<< endl);
408  // Allocate for all threads first
409  m5p_sensitivity = (FLTNB*****)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(FLTNB****));
410  // Allocate the rest only for all threads
411  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
412  {
413  m5p_sensitivity[th] = (FLTNB****)malloc(mp_ID->GetNbTimeFrames()*sizeof(FLTNB***));
414  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
415  {
416  m5p_sensitivity[th][fr] = (FLTNB***)malloc(mp_ID->GetNbRespGates()*sizeof(FLTNB**));
417  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
418  {
419  m5p_sensitivity[th][fr][rg] = (FLTNB**)malloc(mp_ID->GetNbCardGates()*sizeof(FLTNB*));
420  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
421  {
422  m5p_sensitivity[th][fr][rg][cg] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
423  // Initialize to 0.
424  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++) m5p_sensitivity[th][fr][rg][cg][v] = 0.;
425  }
426  }
427  }
428  }
429  // Set the flag for loaded sensitivity to false
430  m_loadedSensitivity = false;
431  }
432 }
433 
434 // =====================================================================
435 // ---------------------------------------------------------------------
436 // ---------------------------------------------------------------------
437 // =====================================================================
438 
440 {
441  // Check the first pointer
442  if (m5p_sensitivity)
443  {
444  // Verbose
445  if (m_verbose>=3) Cout("oImageSpace::DeallocateSensitivityImage() -> Free memory"<< endl);
446  // If loaded sensitivity then fully deallocate only the first thread pointer
448  {
449  int thread_0 = 0;
450  if (m5p_sensitivity[thread_0])
451  {
452  // Loop on time frames
453  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++) if (m5p_sensitivity[thread_0][fr])
454  {
455  // Loop on respiratory gates
456  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++) if (m5p_sensitivity[thread_0][fr][rg])
457  {
458  // Loop on cardiac gates
459  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++) if (m5p_sensitivity[thread_0][fr][rg][cg])
460  {
461  free(m5p_sensitivity[thread_0][fr][rg][cg]);
462  }
463  free(m5p_sensitivity[thread_0][fr][rg]);
464  }
465  free(m5p_sensitivity[thread_0][fr]);
466  }
467  free(m5p_sensitivity[thread_0]);
468  }
469  }
470  // Otherwise, loop on all threads
471  else
472  {
473  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++) if (m5p_sensitivity[th])
474  {
475  // Loop on time frames
476  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++) if (m5p_sensitivity[th][fr])
477  {
478  // Loop on respiratory gates
479  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++) if (m5p_sensitivity[th][fr][rg])
480  {
481  // Loop on cardiac gates
482  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++) if (m5p_sensitivity[th][fr][rg][cg])
483  {
484  free(m5p_sensitivity[th][fr][rg][cg]);
485  }
486  free(m5p_sensitivity[th][fr][rg]);
487  }
488  free(m5p_sensitivity[th][fr]);
489  }
490  free(m5p_sensitivity[th]);
491  }
492  }
493  // Free the first pointer
494  free(m5p_sensitivity);
495  }
496  // And set the pointer to NULL
497  m5p_sensitivity = NULL;
498 }
499 
500 // =====================================================================
501 // ---------------------------------------------------------------------
502 // ---------------------------------------------------------------------
503 // =====================================================================
504 
506 {
507  // In this function, we want to allocate a new miscellaneous image on the double pointer member
508  // m2p_miscellaneousImage, and then we return the pointer to the allocated image.
509 
510  // Update the number of miscellaneous images
512  // Case 1: the current number of image is 0
514  {
515  // Allocate the first pointer
517  }
518  // Case 2: we already have some miscellaneous images
519  else
520  {
521  // Reallocate the first pointer
523  }
524  // Allocate the image
526  // Return the pointer
528 }
529 
530 // =====================================================================
531 // ---------------------------------------------------------------------
532 // ---------------------------------------------------------------------
533 // =====================================================================
534 
536 {
537  // Test the pointer existence
538  if (m2p_miscellaneousImage!=NULL)
539  {
540  // Verbose
541  if (m_verbose>=3) Cout("oImageSpace::DeallocateMiscellaneousImage() -> Free memory"<< endl);
542  // Loop over all miscellaneous images
543  for (int m=0; m<m_nbMiscellaneousImages; m++)
544  {
545  // Test the pointer and free it
547  }
548  // Free the first pointer
550  }
551 }
552 
553 // =====================================================================
554 // ---------------------------------------------------------------------
555 // ---------------------------------------------------------------------
556 // =====================================================================
557 
559 {
560  if(m_verbose>=3) Cout("oImageSpace::InstantiateVisitedVoxelsImage ..."<< endl);
561  mp_visitedVoxelsImage = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
562  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++) mp_visitedVoxelsImage[v] = 0.;
563 }
564 
565 // =====================================================================
566 // ---------------------------------------------------------------------
567 // ---------------------------------------------------------------------
568 // =====================================================================
569 
571 {
572  if(m_verbose>=3) Cout("oImageSpace::InstantiateVisitedVoxelsImage ..."<< endl);
574  mp_visitedVoxelsImage = NULL;
575 }
576 
577 // =====================================================================
578 // ---------------------------------------------------------------------
579 // ---------------------------------------------------------------------
580 // =====================================================================
581 
582 int oImageSpace::InitMultiModalImage(const vector<string>& a_pathToImage)
583 {
584  // Allocate and read only if a file is provided
585  if (!a_pathToImage.empty())
586  {
587  // Allocate memory (if not already done)
588  if (!m_loadedMultiModal)
589  {
590  m2p_multiModalImage = (FLTNB**)malloc(mp_ID->GetNbMultiModalImages()*sizeof(FLTNB*));
591  for (int nb=0; nb<mp_ID->GetNbMultiModalImages(); nb++)
592  {
593  // Verbose
594  if (m_verbose>=3) Cout("oImageSpace::InitMultiModalImage() -> From file '" << a_pathToImage.at(nb) << "'"<< endl);
595  // Open file
596  ifstream image_file(a_pathToImage.at(nb).c_str(), ios::in|ios::binary);
597  if (!image_file.is_open())
598  {
599  Cerr("***** oImageSpace::InitMultiModalImage() -> Input file '" << a_pathToImage.at(nb) << "' is missing or corrupted !" << endl);
600  // Read failure implies that the multiModal image will never be used, so free all the allocated memory
602  return 1;
603  }
604  m2p_multiModalImage[nb] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
605  // Interfile image reading (INTF_LERP_ENABLED = interpolation allowed)
606  if (IntfReadImage(a_pathToImage.at(nb), m2p_multiModalImage[nb], mp_ID, m_verbose, INTF_LERP_ENABLED))
607  {
608  Cerr("***** oImageSpace::InitMultiModalImage() -> Error reading interfile image '" << a_pathToImage.at(nb) << "' !" << endl);
609  // Read failure implies that the multiModal image will never be used, so free all the allocated memory
611  return 1;
612  }
613  // Close file
614  image_file.close();
615  }
616  // Flag as loaded
617  m_loadedMultiModal = true;
618  }
619  }
620  // End
621  return 0;
622 }
623 
624 // =====================================================================
625 // ---------------------------------------------------------------------
626 // ---------------------------------------------------------------------
627 // =====================================================================
628 
630 {
631  // Check the first pointer
633  {
634  // Verbose
635  if (m_verbose>=3) Cout("oImageSpace::DeallocateMultiModalImage() -> Free memory" << endl);
636  for (int nb=0; nb<mp_ID->GetNbMultiModalImages(); nb++) if (m2p_multiModalImage[nb]) free(m2p_multiModalImage[nb]);
637  free(m2p_multiModalImage);
638  }
639  // Reset the pointer to NULL
640  m2p_multiModalImage = NULL;
641 }
642 
643 // =====================================================================
644 // ---------------------------------------------------------------------
645 // ---------------------------------------------------------------------
646 // =====================================================================
647 
648 int oImageSpace::InitMaskImage(const string& a_pathToImage)
649 {
650  // Allocate and read only if a file is provided
651  if (!a_pathToImage.empty())
652  {
653  // Verbose
654  if (m_verbose>=3) Cout("oImageSpace::InitMaskImage() -> From file '" << a_pathToImage << "'"<< endl);
655  // Open file
656  ifstream image_file(a_pathToImage.c_str(), ios::in|ios::binary);
657  if (!image_file.is_open())
658  {
659  Cerr("***** oImageSpace::InitMaskImage() -> Input file '" << a_pathToImage << "' is missing or corrupted !" << endl);
660  return 1;
661  }
662  // Allocate memory (if not already done)
663  if (!m_loadedMask)
664  {
665  mp_maskImage = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
666  }
667  // Interfile image reading (INTF_LERP_ENABLED = interpolation allowed)
669  {
670  Cerr("***** oImageSpace::InitMaskImage() -> Error reading interfile image '" << a_pathToImage << "' !" << endl);
671  // Read failure implies that the mask image will never be used, so free all the allocated memory
673  return 1;
674  }
675  // Flag as loaded
676  m_loadedMask = true;
677  // Close file
678  image_file.close();
679  }
680  // End
681  return 0;
682 }
683 
684 // =====================================================================
685 // ---------------------------------------------------------------------
686 // ---------------------------------------------------------------------
687 // =====================================================================
688 
690 {
691  // Check pointer
692  if (mp_maskImage)
693  {
694  // Verbose
695  if (m_verbose>=3) Cout("oImageSpace::DeallocateMaskImage() -> Free memory" << endl);
696  free(mp_maskImage);
697  }
698  // Set the pointer to NULL
699  mp_maskImage = NULL;
700  // flag as not loaded
701  m_loadedMask = false;
702 }
703 
704 
705 
706 
707 // =====================================================================
708 // ---------------------------------------------------------------------
709 // ---------------------------------------------------------------------
710 // =====================================================================
711 /*
712  \fn InstantiateOutputImage()
713  \brief Instanciate Image matrices dedicated to output writing on disk
714  \details Additionnal output image matrix is needed in the case the reconstruction uses intrinsic temporal basis functions
715  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.
716  If no intrinsic temporal basis functions are used, m4p_outputImage just refers to the address of the image matrix containing the regular image.
717 */
719 {
720  if(m_verbose>=3) Cout("oImageSpace::InstantiateOutputImage ..."<< endl);
721 
722  // Here, we will allocate the first 3 pointers to the dynamic dimensions
724  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
725  {
726  m4p_outputImage[fr] = new FLTNB**[mp_ID->GetNbRespGates()];
727  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
728  {
729  m4p_outputImage[fr][rg] = new FLTNB*[mp_ID->GetNbCardGates()];
730  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
731  {
732  // For the last pointer (i.e. the voxels), we allocate them only if we strictly need them, which
733  // means in the case where one of the dynamic dimensions makes internal use of basis functions,
734  // otherwise, we will just make this output image pointing to the forward image which is useless
735  // at the moment of computing the output image.
736  if ( mp_ID->GetTimeStaticFlag() &&
739  {
740  // In this case, the time frames and basis functions represent the same thing, same for respiratory/cardiac
741  // gates and their equivalent basis functions
742  m4p_outputImage[fr][rg][cg] = m4p_forwardImage[fr][rg][cg];
743  }
744  else
745  {
746  // We allocate the output image
747  m4p_outputImage[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
748  }
749  }
750  }
751  }
752 }
753 
754 
755 
756 
757 // =====================================================================
758 // ---------------------------------------------------------------------
759 // ---------------------------------------------------------------------
760 // =====================================================================
761 /*
762  \fn DeallocateOutputImage()
763  \brief Free memory for the Image matrices dedicated to output writing on disk
764 */
766 {
767  if(m_verbose>=3) Cout("oImageSpace::DeallocateOutputImage ..."<< endl);
768 
769  if (m4p_outputImage)
770  {
771  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
772  {
773  if (m4p_outputImage[fr])
774  {
775  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
776  {
777  if (m4p_outputImage[fr][rg])
778  {
779  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
780  {
781  // Distinguish between intrinsic use of dynamic basis functions or not (in the latter,
782  // the m4p_outputImage points to the m4p_forwardImage, so we do not touch it)
783  if ( (!mp_ID->GetTimeStaticFlag() ||
784  !mp_ID->GetRespStaticFlag() ||
785  !mp_ID->GetCardStaticFlag()) &&
786  m4p_outputImage[fr][rg][cg] )
787  {
788  delete m4p_outputImage[fr][rg][cg];
789  }
790  }
791  delete[] m4p_outputImage[fr][rg];
792  }
793  }
794  delete[] m4p_outputImage[fr];
795  }
796  }
797  delete[] m4p_outputImage;
798  }
799  m4p_outputImage = NULL;
800 }
801 
802 
803 
804 // =====================================================================
805 // ---------------------------------------------------------------------
806 // ---------------------------------------------------------------------
807 // =====================================================================
808 /*
809  \fn InstantiateBwdImageForDeformation()
810  \brief Memory allocation for the buffer backward image required for image-based deformation
811 
812 void oImageSpace::InstantiateRefImagesForDeformation()
813 {
814  if(m_verbose>=3) Cout("oImageSpace::InstantiateBwdImageForDeformation ..."<< endl);
815 
816  m5p_refDynBackwardImage = new FLTNB****[m_nbBackwardImages];
817 
818  for (int img=0; img<m_nbBackwardImages; img++)
819  {
820  m5p_refDynBackwardImage[img] = new FLTNB***[mp_ID->GetNbTimeBasisFunctions()];
821  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
822  {
823  // The fourth pointer is for the number of respiratory basis functions
824  m5p_refDynBackwardImage[img][tbf] = new FLTNB**[mp_ID->GetNbRespBasisFunctions()];
825  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
826  {
827  // The fifth pointer is for the number of cardiac basis functions
828  m5p_refDynBackwardImage[img][tbf][rbf] = new FLTNB*[mp_ID->GetNbCardBasisFunctions()];
829  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
830  {
831  // The sixth pointer is for the 3D space
832  m5p_refDynBackwardImage[img][tbf][rbf][cbf] = new FLTNB[mp_ID->GetNbVoxXYZ()];
833  }
834  }
835  }
836  }
837 
838 }
839 */
840 
841 /*
842  \fn oImageSpace::InstantiateRefImagesForDeformation
843  \brief Allocate memory for the buffer sensitivity image required for image-based deformation. This function is called from the Deformation Manager
844 */
846 {
847  if(m_verbose>=3) Cout("oImageSpace::InstantiateRefImagesForDeformation ..."<< endl);
848 
849 
850  // Forward image
852  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
853  {
854  // The fourth pointer is for the number of respiratory basis functions
856  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
857  {
858  // The fifth pointer is for the number of cardiac basis functions
860  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
861  {
862  // The sixth pointer is for the 3D space
863  m4p_refDynForwardImage[tbf][rbf][cbf] = new FLTNB[mp_ID->GetNbVoxXYZ()];
864  }
865  }
866  }
867 
868 
869  // Backward image
871 
872  for (int img=0; img<m_nbBackwardImages; img++)
873  {
875  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
876  {
877  // The fourth pointer is for the number of respiratory basis functions
879  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
880  {
881  // The fifth pointer is for the number of cardiac basis functions
883  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
884  {
885  // The sixth pointer is for the 3D space
886  m5p_refDynBackwardImage[img][tbf][rbf][cbf] = new FLTNB[mp_ID->GetNbVoxXYZ()];
887  }
888  }
889  }
890  }
891 
892 
893  // Sensitivity image
894  if(m_loadedSensitivity == false)
895  {
897 
898  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
899  {
901  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
902  {
904  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
905  m4p_refDynSensitivityImage[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
906  }
907  }
908  }
909 
910 }
911 
912 
913 // =====================================================================
914 // ---------------------------------------------------------------------
915 // ---------------------------------------------------------------------
916 // =====================================================================
917 /*
918  \fn DeallocateBwdImageForDeformation()
919  \brief Free memory for the buffer backward image required for image-based deformation
920 
921 void oImageSpace::DeallocateBwdImageForDeformation()
922 {
923  if(m_verbose>=3) Cout("oImageSpace::DeallocateBwdImageForDeformation ..."<< endl);
924 
925  for (int img=0; img<m_nbBackwardImages; img++)
926  {
927  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
928  {
929  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
930  {
931  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
932  {
933  if (m5p_refDynBackwardImage[img][tbf][rbf][cbf]) delete m5p_refDynBackwardImage[img][tbf][rbf][cbf];
934  }
935  if (m5p_refDynBackwardImage[img][tbf][rbf]) delete m5p_refDynBackwardImage[img][tbf][rbf];
936  }
937  if (m5p_refDynBackwardImage[img][tbf]) delete[] m5p_refDynBackwardImage[img][tbf];
938  }
939  if (m5p_refDynBackwardImage[img]) delete[] m5p_refDynBackwardImage[img];
940  }
941  if (m5p_refDynBackwardImage) delete[] m5p_refDynBackwardImage;
942  m5p_refDynBackwardImage = NULL;
943 }
944 */
945 
946 
947 
948 
949 // =====================================================================
950 // ---------------------------------------------------------------------
951 // ---------------------------------------------------------------------
952 // =====================================================================
953 /*
954  \fn oImageSpace::DeallocateRefImagesForDeformation
955  \brief Free memory for the buffer sensitivity image required for image-based deformation. This function is called from the Deformation Manager
956 */
958 {
959  if(m_verbose>=3) Cout("oImageSpace::DeallocateRefImagesForDeformation ..."<< endl);
960 
961  // Forward image
962  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
963  {
964  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
965  {
966  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
967  {
968  if (m4p_refDynForwardImage[tbf][rbf][cbf]) delete[] m4p_refDynForwardImage[tbf][rbf][cbf];
969  }
970  if (m4p_refDynForwardImage[tbf][rbf]) delete[] m4p_refDynForwardImage[tbf][rbf];
971  }
972  if (m4p_refDynForwardImage[tbf]) delete[] m4p_refDynForwardImage[tbf];
973  }
975 
976  m4p_refDynForwardImage = NULL;
977 
978  // Backward image
979  for (int img=0; img<m_nbBackwardImages; img++)
980  {
981  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
982  {
983  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
984  {
985  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
986  {
987  if (m5p_refDynBackwardImage[img][tbf][rbf][cbf]) delete[] m5p_refDynBackwardImage[img][tbf][rbf][cbf];
988  }
989  if (m5p_refDynBackwardImage[img][tbf][rbf]) delete[] m5p_refDynBackwardImage[img][tbf][rbf];
990  }
991  if (m5p_refDynBackwardImage[img][tbf]) delete[] m5p_refDynBackwardImage[img][tbf];
992  }
993  if (m5p_refDynBackwardImage[img]) delete[] m5p_refDynBackwardImage[img];
994  }
997 
998  // Sensitivity image
999  if(m_loadedSensitivity == false)
1000  {
1001  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1002  {
1003  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1004  {
1005  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1006  if (m4p_refDynSensitivityImage[fr][rg][cg]) delete[] m4p_refDynSensitivityImage[fr][rg][cg];
1007 
1008  if (m4p_refDynSensitivityImage[fr][rg]) delete[] m4p_refDynSensitivityImage[fr][rg];
1009  }
1011  }
1014  }
1015 }
1016 
1017 // =====================================================================
1018 // ---------------------------------------------------------------------
1019 // ---------------------------------------------------------------------
1020 // =====================================================================
1021 /*
1022  \fn InstantiateSensImageForDeformation()
1023  \brief Memory allocation for the buffer sensitivity image required for image-based deformation (only for PET HISTOGRAM mode)
1024 
1025 void oImageSpace::InstantiateSensImageForDeformation()
1026 {
1027  if(m_verbose>=3) Cout("oImageSpace::InstantiateSensImageForDeformation ..."<< endl);
1028 
1029  m4p_refDynSensitivityImage = new FLTNB***[mp_ID->GetNbTimeFrames()];
1030 
1031  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1032  {
1033  m4p_refDynSensitivityImage[fr] = new FLTNB**[mp_ID->GetNbRespGates()];
1034  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1035  {
1036  m4p_refDynSensitivityImage[fr][rg] = new FLTNB*[mp_ID->GetNbCardGates()];
1037  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1038  m4p_refDynSensitivityImage[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
1039  }
1040  }
1041 }
1042 */
1043 
1044 
1045 
1046 // =====================================================================
1047 // ---------------------------------------------------------------------
1048 // ---------------------------------------------------------------------
1049 // =====================================================================
1050 /*
1051  \fn DeallocateBwdImageForDeformation()
1052  \brief Free memory for the buffer sensitivity image required for image-based deformation
1053 
1054 void oImageSpace::DeallocateSensImageForDeformation()
1055 {
1056  if(m_verbose>=3) Cout("oImageSpace::DeallocateSensImageForDeformation ..."<< endl);
1057 
1058  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1059  {
1060  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1061  {
1062  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1063  if (m4p_refDynSensitivityImage[fr][rg][cg]) delete[] m4p_refDynSensitivityImage[fr][rg][cg];
1064 
1065  if (m4p_refDynSensitivityImage[fr][rg]) delete[] m4p_refDynSensitivityImage[fr][rg];
1066  }
1067  if (m4p_refDynSensitivityImage[fr]) delete[] m4p_refDynSensitivityImage[fr];
1068  }
1069  if (m4p_refDynSensitivityImage) delete[] m4p_refDynSensitivityImage;
1070  m4p_refDynSensitivityImage = NULL;
1071 }
1072 */
1073 
1074 
1075 // =====================================================================
1076 // ---------------------------------------------------------------------
1077 // ---------------------------------------------------------------------
1078 // =====================================================================
1079 /*
1080  \fn LMS_DeallocateAttenuationImage()
1081  \brief Free memory for the Attenuation image matrices (for analytical projection or list-mode sensitivity generation)
1082 */
1084 {
1085  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateAttenuationImage ..."<< endl);
1086 
1087  if(m4p_attenuation != NULL)
1088  {
1089  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1090  {
1091  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
1092  {
1093  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
1094  if(m4p_attenuation[fr][rg][cg] != NULL) delete m4p_attenuation[fr][rg][cg];
1095 
1096  if(m4p_attenuation[fr][rg] != NULL) delete m4p_attenuation[fr][rg];
1097  }
1098  if(m4p_attenuation[fr] != NULL) delete m4p_attenuation[fr];
1099  }
1100  if(m4p_attenuation != NULL) delete[] m4p_attenuation;
1101  m4p_attenuation = NULL;
1102  }
1103 }
1104 
1105 
1106 
1107 
1108 // =====================================================================
1109 // ---------------------------------------------------------------------
1110 // ---------------------------------------------------------------------
1111 // =====================================================================
1112 
1113 int oImageSpace::InitAttenuationImage(const string& a_pathToAtnImage)
1114 {
1115  // todo : read the number of dynamic images from Interfile header
1116  // Allocate memory and initialize matrices according to the number of fr/rg/cg
1117 
1118  // Allocate only if an image has been provided
1119  if (a_pathToAtnImage.empty()) return 0;
1120  // Verbose
1121  if (m_verbose>=VERBOSE_NORMAL) Cout("oImageSpace::InitAttenuationImage() -> Initialize and read attenuation image from file '" << a_pathToAtnImage << "'"<< endl);
1122  // Allocate only if not already done
1123  if (!m4p_attenuation)
1124  {
1126  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1127  {
1129  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
1130  {
1131  m4p_attenuation[fr][rg] = new FLTNB*[mp_ID->GetNb2ndMotImgsForLMS()];
1132 
1133  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
1134  m4p_attenuation[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
1135  }
1136  }
1137  }
1138 
1139  // Open file
1140  ifstream image_file(a_pathToAtnImage.c_str(), ios::in|ios::binary); // Read the corresponding attenuation image
1141  // Check opening
1142  if (!image_file.is_open())
1143  {
1144  Cerr("***** oImageSpace::InitAttenuationImage() -> Failed to open attenuation image from file '" << a_pathToAtnImage << "' !" << endl);
1145  return 1;
1146  }
1147  else
1148  {
1149  // Interfile image reading (INTF_LERP_ENABLED = interpolation allowed)
1150  if (IntfReadImage(a_pathToAtnImage, m4p_attenuation, mp_ID, m_verbose, INTF_LERP_ENABLED) )
1151  {
1152  Cerr("***** oImageSpace::InitAttenuationImage() -> An error occurred while reading from file '" << a_pathToAtnImage << "' !" << endl);
1153  return 1;
1154  }
1155 
1156  // Check if we have to copy the attenuation images to potential other frames/gates
1157 
1158  // Recover the number of dynamic images in the attenuation imag file
1159  int nb_atn_frames = 1, nb_atn_rgates = 1, nb_atn_cgates = 1;
1160 
1161  if( (IntfKeyGetValueFromFile(a_pathToAtnImage, "number of time frames := " , &nb_atn_frames, 1, KEYWORD_OPTIONAL) == KEYWORD_OPTIONAL_ERROR )
1162  ||(IntfKeyGetValueFromFile(a_pathToAtnImage, "number of respiratory gates := ", &nb_atn_rgates, 1, KEYWORD_OPTIONAL) == KEYWORD_OPTIONAL_ERROR )
1163  ||(IntfKeyGetValueFromFile(a_pathToAtnImage, "number of cardiac gates := " , &nb_atn_cgates, 1, KEYWORD_OPTIONAL) == KEYWORD_OPTIONAL_ERROR ) )
1164  {
1165  Cerr("***** oImageSpace::InitAttenuationImage() -> An error occurred while reading'number of time frames', 'number of respiratory gates' or 'number of cardiac gates' keys from file '" << a_pathToAtnImage << "' !" << endl);
1166  return 1;
1167  }
1168 
1169 
1170  // If more than one frame, Copy attenuation image to each other frames
1171 
1172  if(nb_atn_cgates == 1 && mp_ID->GetNb2ndMotImgsForLMS()>1)
1173  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1174  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
1175  for(int cg=1 ; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
1176  for(int v=0 ; v<mp_ID->GetNbVoxXYZ(); v++)
1177  m4p_attenuation[fr][rg][cg][v] = m4p_attenuation[fr][rg][0][v];
1178 
1179 
1180 
1181  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1182  if(nb_atn_rgates == 1 && mp_ID->GetNb1stMotImgsForLMS(fr)>1)
1183  for(int rg=1 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
1184  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
1185  for(int v=0 ; v<mp_ID->GetNbVoxXYZ(); v++)
1186  m4p_attenuation[fr][rg][cg][v] = m4p_attenuation[fr][0][cg][v];
1187 
1188 
1189 
1190  if(nb_atn_frames == 1 && mp_ID->GetNbTimeFrames()>1)
1191  for(int fr=1 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1192  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
1193  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
1194  for(int v=0 ; v<mp_ID->GetNbVoxXYZ(); v++)
1195  // If chronological motion correction is enabled, the number of motion image may be different for each frame
1196  // Check that here, just copy the attenuation image of the first motion image if it is the case
1197  m4p_attenuation[fr][rg][cg][v] = ( mp_ID->GetNb1stMotImgsForLMS(fr) == mp_ID->GetNb1stMotImgsForLMS(0) ) ?
1198  m4p_attenuation[0][rg][cg][v] :
1199  m4p_attenuation[0][0][cg][v] ;
1200 
1201 
1202  }
1203 
1204  // Close file
1205  image_file.close();
1206  return 0;
1207 }
1208 
1209 // =====================================================================
1210 // ---------------------------------------------------------------------
1211 // ---------------------------------------------------------------------
1212 // =====================================================================
1213 /*
1214  \fn InitImage
1215  \param a_pathToInitialImage : path to an existing image
1216  \param a_value : value to initialize each voxel with, if an input image is not provided
1217  \brief Initialize the main image, either using:
1218  - an existing image at path 'a_pathToInitialImage'
1219  - initialize each voxel with 'a_value'.
1220  \return 0 if success, positive value otherwise
1221 */
1222 int oImageSpace::InitImage(const string& a_pathToInitialImage, FLTNB a_value)
1223 {
1224  if(m_verbose>=3) Cout("oImageSpace::InitImage ..."<< endl);
1225 
1226  if (!a_pathToInitialImage.empty()) // Image initiale has been provided
1227  {
1228  if (LoadInitialImage(a_pathToInitialImage) )
1229  {
1230  Cerr("***** oImageSpace::InitImage()-> Error while trying to read file at " << a_pathToInitialImage << endl);
1231  return 1;
1232  }
1233  }
1234 
1235  else // Uniform initialization
1236  {
1237  // We initialize the content of the cylindrical FOV with the provided
1238  // value, and the exterior with a ratio of the provided value
1239  FLTNB exterior_fov_value = a_value * 1.; // For the moment we let the same value everywhere
1240  // Precast half the number of voxels over X and Y minus 1 (for efficiency)
1241  FLTNB flt_base_x = 0.5*((FLTNB)(mp_ID->GetNbVoxX()-1));
1242  FLTNB flt_base_y = 0.5*((FLTNB)(mp_ID->GetNbVoxY()-1));
1243  // Compute FOV elipse radius over X and Y, then squared
1244  FLTNB squared_radius_x = 0.5 * ((FLTNB)(mp_ID->GetNbVoxX())) * mp_ID->GetVoxSizeX();
1245  squared_radius_x *= squared_radius_x;
1246  FLTNB squared_radius_y = 0.5 * ((FLTNB)(mp_ID->GetNbVoxY())) * mp_ID->GetVoxSizeY();
1247  squared_radius_y *= squared_radius_y;
1248  // We assume that the computation of the distance from the center for a given
1249  // voxel and comparing it with the output FOV percentage costs more than performing
1250  // the loops in an inverse order compared to how the image is stored in memory.
1251  // Thus we begin the loops over X, then Y, then we test and if test passes, we
1252  // do the remaining loop over Z and over all dynamic dimensions.
1253  int x;
1254  #pragma omp parallel for private(x) schedule(guided)
1255  for (x=0; x<mp_ID->GetNbVoxX(); x++)
1256  {
1257  // Compute X distance from image center, then squared
1258  FLTNB squared_distance_x = (((FLTNB)x)-flt_base_x) * mp_ID->GetVoxSizeX();
1259  squared_distance_x *= squared_distance_x;
1260  // Loop over Y
1261  for (int y=0; y<mp_ID->GetNbVoxY(); y++)
1262  {
1263  // Compute Y distance from image center, then squared
1264  FLTNB squared_distance_y = (((FLTNB)y)-flt_base_y) * mp_ID->GetVoxSizeY();
1265  squared_distance_y *= squared_distance_y;
1266  // The value to affect to the voxel
1267  FLTNB affected_value = 0.;
1268  // Test if the voxel is inside the FOV elipse, then we affect the provided value
1269  if ( squared_distance_x/squared_radius_x + squared_distance_y/squared_radius_y <= 1. ) affected_value = a_value;
1270  // Else if outside, we affect the exterior value
1271  else affected_value = exterior_fov_value;
1272  // Loop over Z
1273  for (int z=0; z<mp_ID->GetNbVoxZ(); z++)
1274  {
1275  // Compute global voxel index
1276  INTNB index = z*mp_ID->GetNbVoxXY() + y*mp_ID->GetNbVoxX() + x;
1277  // Loops on dynamic dimensions
1278  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1279  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1280  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1281  // Affect image value
1282  m4p_image[tbf][rbf][cbf][index] = affected_value;
1283  }
1284  }
1285  }
1286  }
1287 
1288  return 0;
1289 }
1290 
1291 
1292 
1293 
1294 // =====================================================================
1295 // ---------------------------------------------------------------------
1296 // ---------------------------------------------------------------------
1297 // =====================================================================
1298 /*
1299  \fn LoadInitialImage()
1300  \param a_pathToImage : path to an existing image
1301  \brief Load the initial image provided by the user in the corresponding matrix
1302  \return 0 if success, positive value otherwise
1303 */
1304 int oImageSpace::LoadInitialImage(const string& a_pathToImage)
1305 {
1306  if(m_verbose>=3) Cout("oImageSpace::LoadInitialImage ..."<< endl);
1307 
1308  ifstream image_file;
1309  image_file.open(a_pathToImage.c_str(), ios::in | ios::binary);
1310 
1311  if(!image_file.is_open())
1312  {
1313  Cerr("***** oImageSpace::LoadInitialImage()-> Error reading file!" << endl);
1314  return 1;
1315  }
1316  else
1317  {
1318  // Interfile image reading (INTF_LERP_DISABLED = no interpolation allowed)
1319  // if(IntfReadImage(a_pathToImage, m4p_image, mp_ID, m_verbose, INTF_LERP_DISABLED))
1321  {
1322  Cerr("***** oImageSpace::LoadInitialImage()-> Error reading Interfile : " << a_pathToImage << " !" << endl);
1323  return 1;
1324  }
1325  }
1326  image_file.close();
1327 
1328  return 0;
1329 }
1330 
1331 
1332 
1333 
1334 // =====================================================================
1335 // ---------------------------------------------------------------------
1336 // ---------------------------------------------------------------------
1337 // =====================================================================
1338 /*
1339  \fn InitializationBackwardImage()
1340  \brief Initialize each voxel of the backward images to 0, also for sensitivity if not loaded (estimated on the fly).
1341 */
1343 {
1344  if(m_verbose>=3) Cout("oImageSpace::InitBackwardImage ..." << endl);
1345 
1346  // Reset backward images to 0.
1347  for (int img=0; img<m_nbBackwardImages; img++)
1348  {
1349  int th;
1350  #pragma omp parallel for private(th) schedule(static, 1)
1351  for (th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
1352  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1353  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1354  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1355  for (INTNB v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1356  {
1357  m6p_backwardImage[img][th][tbf][rbf][cbf][v] = 0.;
1358  }
1359  }
1360 
1361  // If on-the-fly sensitivity, then reset to 0.
1362  if (!m_loadedSensitivity)
1363  {
1364  int th;
1365  #pragma omp parallel for private(th) schedule(static, 1)
1366  for (th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
1367  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1368  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1369  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1370  for (INTNB v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1371  m5p_sensitivity[th][fr][rg][cg][v] = 0.;
1372  }
1373 }
1374 
1375 
1376 
1377 
1378 
1379 
1380 // =====================================================================
1381 // ---------------------------------------------------------------------
1382 // ---------------------------------------------------------------------
1383 // =====================================================================
1384 /*
1385  \fn InitSensitivityImage()
1386  \param a_pathToSensitivityImage : path to the sensitivity image (should be provided only in list-mode reconstruction)
1387  \brief Initialization for the sensitivity image matrices
1388  \details Sensitivity image initialization depends on the reconstruction mode :
1389  - list-mode: Sensitivity image has been computed before reconstruction, it is loaded from the path provided in parameter
1390  - histogram: Sensitivity image is calculated on the fly during reconstruction.
1391  First dimension (thread) is only used in histogram mode, as the on-the-fly sensitivity image computation must be thread safe
1392  \return 0 if success, positive value otherwise
1393 */
1394 int oImageSpace::InitSensitivityImage(const string& a_pathToSensitivityImage)
1395 {
1396  if(m_verbose>=3) Cout("oImageSpace::InitSensitivityImage ..."<< endl);
1397 
1398  // =====================================================================================================
1399  // Case 1: a path to a sensitivity image is provided, meaning that we are in listmode and do not compute
1400  // the sensitivity on-the-fly. In that case, we do not take the multi-threading into account.
1401  // =====================================================================================================
1402 
1403  if (!a_pathToSensitivityImage.empty())
1404  {
1405  // Open file
1406  ifstream input_file;
1407  input_file.open(a_pathToSensitivityImage.c_str(), ios::binary | ios::in);
1408  // Read file
1409  if (input_file.is_open())
1410  {
1411  // Interfile image reading (INTF_LERP_DISABLED = no interpolation allowed)
1412  if(IntfReadImage(a_pathToSensitivityImage, m5p_sensitivity[0], mp_ID, m_verbose, INTF_LERP_DISABLED) )
1413  {
1414  Cerr("***** oImageSpace::InitSensitivityImage()-> Error reading Interfile : " << a_pathToSensitivityImage << " !" << endl);
1415  return 1;
1416  }
1417  input_file.close();
1418  m_loadedSensitivity = true;
1419  }
1420  else
1421  {
1422  Cerr("***** oImageSpace::InitSensitivityImage() -> Input sensitivity file '" << a_pathToSensitivityImage << "' is missing or corrupted !" << endl);
1423  return 1;
1424  }
1425  }
1426 
1427  // =====================================================================================================
1428  // Case 2: no sensitivity provided, thus we are in histogram mode and will commpute it on-the-fly. We
1429  // thus need it to be thread-safe, so we allocated for all threads.
1430  // =====================================================================================================
1431 
1432  else
1433  {
1434 /*
1435  // Standard initialization
1436  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
1437  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1438  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1439  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1440  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1441  m5p_sensitivity[th][fr][rg][cg][v] = -1.;
1442 * */
1443  // No need to initialize to any value, as it will be done by the InitBackwardImage function which is called at the beginning
1444  // of each subset, and deals with the sensitivity initialization to 0 when m_loadedSensitivity is false
1445  m_loadedSensitivity = false;
1446  }
1447 
1448  // End
1449  return 0;
1450 }
1451 
1452 
1453 
1454 
1455 
1456 // =====================================================================
1457 // ---------------------------------------------------------------------
1458 // ---------------------------------------------------------------------
1459 // =====================================================================
1460 /*
1461  \fn oImageSpace::InitRefImagesForDeformation
1462  \brief Initialize the references image dedicated to image-based deformation. This function is called from the Deformation Manager
1463 */
1465 {
1466  if(m_verbose>=3) Cout("oDeformationManager::InitBwdImageForDeformation ..." << endl);
1467 
1468  // forward image
1469  for (int img=0; img<m_nbBackwardImages; img++)
1470  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1471  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1472  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1473  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1474  m4p_refDynForwardImage[tbf][rbf][cbf][v] = m4p_forwardImage[tbf][rbf][cbf][v];
1475 
1476  // backward image
1477  for (int img=0; img<m_nbBackwardImages; img++)
1478  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1479  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1480  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1481  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1482  m5p_refDynBackwardImage[img][tbf][rbf][cbf][v] = 0.;
1483 
1484  if(m_loadedSensitivity == false)
1485  {
1486  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1487  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1488  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1489  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
1490  m4p_refDynSensitivityImage[fr][rg][cg][v] = 0.;
1491  }
1492 }
1493 
1494 
1495 
1496 
1497 
1498 // =====================================================================
1499 // ---------------------------------------------------------------------
1500 // ---------------------------------------------------------------------
1501 // =====================================================================
1502 /*
1503  \fn InitializationSensImageForDeformation()
1504  \brief Initialize the buffer sensitivity image dedicated to image-based deformation, if required (histogram mode, as sensitivity is not loaded)
1505 
1506 void oImageSpace::InitSensImageForDeformation()
1507 {
1508  if(m_verbose>=3) Cout("oImageSpace::InitSensImageForDeformation ..."<< endl);
1509 
1510  if(m_loadedSensitivity == false)
1511  {
1512  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1513  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1514  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1515  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
1516  m4p_refDynSensitivityImage[fr][rg][cg][v] = 0.;
1517  }
1518 }
1519 */
1520 
1521 
1522 
1523 
1524 // =====================================================================
1525 // ---------------------------------------------------------------------
1526 // ---------------------------------------------------------------------
1527 // =====================================================================
1528 /*
1529  \fn InstantiateOutputImage()
1530  \brief Compute output image using the m4p_image matrix and the time/respiratory/cardiac basis functions. Store the result in the m4p_outputImage matrix
1531  \details If time/respiratory/cardiac basis functions have been initialized, this function has no effect.
1532 */
1534 {
1535  if(m_verbose>=3) Cout("oImageSpace::ComputeOutputImage ..."<< endl);
1536 
1537  // First loop on frames
1538  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1539  {
1540  // Second loop on respiratory gates
1541  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1542  {
1543  // Third loop on cardiac gates
1544  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1545  {
1546  // Reset m4p_outputImage to 0
1547  int v;
1548  #pragma omp parallel for private(v) schedule(guided)
1549  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1550  m4p_outputImage[fr][rg][cg][v] = 0.;
1551  // First loop on time basis functions
1552  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1553  {
1554  // Get frame/basis coefficient and continue if null
1555  FLTNB time_basis_coef = mp_ID->GetTimeBasisCoefficient(tbf,fr);
1556  if (time_basis_coef==0.) continue;
1557  // Second loop on respiratory basis functions
1558  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1559  {
1560  // Get resp_gate/basis coefficient and continue if null
1561  FLTNB resp_basis_coef = mp_ID->GetRespBasisCoefficient(rbf,rg);
1562  if (resp_basis_coef==0.) continue;
1563  // Third loop on cardiac basis functions
1564  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1565  {
1566  // Get card_gate_basis coefficient and continue if null
1567  FLTNB card_basis_coef = mp_ID->GetCardBasisCoefficient(cbf,cg);
1568  if (card_basis_coef==0.) continue;
1569  // Compute global coefficient
1570  FLTNB global_basis_coeff = time_basis_coef * resp_basis_coef * card_basis_coef;
1571  // Loop on voxel with OpenMP (let the chunk size by default as all values are aligned in memory)
1572  #pragma omp parallel for private(v) schedule(guided)
1573  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1574  {
1575  // Add contribution from these basis functions
1576  m4p_outputImage[fr][rg][cg][v] += m4p_image[tbf][rbf][cbf][v] * global_basis_coeff;
1577  }
1578  }
1579  }
1580  }
1581  }
1582  }
1583  }
1584 }
1585 
1586 // =====================================================================
1587 // ---------------------------------------------------------------------
1588 // ---------------------------------------------------------------------
1589 // =====================================================================
1590 
1592 {
1593  // Note: The output image is copied in the forward image at this step
1594 
1595  // If no flip, then return
1596  if ( !mp_ID->GetFlipOutX() &&
1597  !mp_ID->GetFlipOutY() &&
1598  !mp_ID->GetFlipOutZ() ) return 0;
1599  // Verbose
1600  if (m_verbose>=2)
1601  {
1602  Cout("oImageSpace::ApplyOutputFlip() -> Flip image" << endl);
1603  if (mp_ID->GetFlipOutX()) Cout(" --> Over X" << endl);
1604  if (mp_ID->GetFlipOutY()) Cout(" --> Over Y" << endl);
1605  if (mp_ID->GetFlipOutZ()) Cout(" --> Over Z" << endl);
1606  }
1607 /*
1608  // First loop on frames
1609  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1610  {
1611  // Second loop on respiratory gates
1612  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1613  {
1614  // Third loop on cardiac gates
1615  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1616  {
1617 */
1618  // First loop on time basis functions
1619  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1620  {
1621  // Second loop on respiratory basis functions
1622  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1623  {
1624  // Third loop on cardiac basis functions
1625  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1626  {
1627  // Flip over Z
1628  if (mp_ID->GetFlipOutZ())
1629  {
1630  // Half loop over Z
1631  for (INTNB z_1=0; z_1<mp_ID->GetNbVoxZ()/2; z_1++)
1632  {
1633  // Compute opposite Z
1634  INTNB z_2 = mp_ID->GetNbVoxZ() - 1 - z_1;
1635  // For efficiency
1636  INTNB base_z1 = z_1 * mp_ID->GetNbVoxXY();
1637  INTNB base_z2 = z_2 * mp_ID->GetNbVoxXY();
1638  // Loop over Y
1639  for (INTNB y=0; y<mp_ID->GetNbVoxY(); y++)
1640  {
1641  // For efficiency
1642  INTNB base_y = y * mp_ID->GetNbVoxX();
1643  // Loop over X
1644  for (INTNB x=0; x<mp_ID->GetNbVoxX(); x++)
1645  {
1646  // Compute both indices
1647  INTNB indice_1 = base_z1 + base_y + x;
1648  INTNB indice_2 = base_z2 + base_y + x;
1649  // Switch voxels
1650 /*
1651  FLTNB buffer = m4p_outputImage[fr][rg][cg][indice_1];
1652  m4p_outputImage[fr][rg][cg][indice_1] = m4p_outputImage[fr][rg][cg][indice_2];
1653  m4p_outputImage[fr][rg][cg][indice_2] = buffer;
1654 */
1655  FLTNB buffer = m4p_forwardImage[tbf][rbf][cbf][indice_1];
1656  m4p_forwardImage[tbf][rbf][cbf][indice_1] = m4p_forwardImage[tbf][rbf][cbf][indice_2];
1657  m4p_forwardImage[tbf][rbf][cbf][indice_2] = buffer;
1658  }
1659  }
1660  }
1661  }
1662  // Flip over Y
1663  if (mp_ID->GetFlipOutY())
1664  {
1665  // Loop over Z
1666  for (INTNB z=0; z<mp_ID->GetNbVoxZ(); z++)
1667  {
1668  // For efficiency
1669  INTNB base_z = z * mp_ID->GetNbVoxXY();
1670  // Half loop over Y
1671  for (INTNB y_1=0; y_1<mp_ID->GetNbVoxY()/2; y_1++)
1672  {
1673  // Compute opposite Y
1674  INTNB y_2 = mp_ID->GetNbVoxY() - 1 - y_1;
1675  // For efficiency
1676  INTNB base_y1 = y_1 * mp_ID->GetNbVoxX();
1677  INTNB base_y2 = y_2 * mp_ID->GetNbVoxX();
1678  // Loop over X
1679  for (INTNB x=0; x<mp_ID->GetNbVoxX(); x++)
1680  {
1681  // Compute both indices
1682  INTNB indice_1 = base_z + base_y1 + x;
1683  INTNB indice_2 = base_z + base_y2 + x;
1684  // Switch voxels
1685 /*
1686  FLTNB buffer = m4p_outputImage[fr][rg][cg][indice_1];
1687  m4p_outputImage[fr][rg][cg][indice_1] = m4p_outputImage[fr][rg][cg][indice_2];
1688  m4p_outputImage[fr][rg][cg][indice_2] = buffer;
1689 */
1690  FLTNB buffer = m4p_forwardImage[tbf][rbf][cbf][indice_1];
1691  m4p_forwardImage[tbf][rbf][cbf][indice_1] = m4p_forwardImage[tbf][rbf][cbf][indice_2];
1692  m4p_forwardImage[tbf][rbf][cbf][indice_2] = buffer;
1693  }
1694  }
1695  }
1696  }
1697  // Flip over X
1698  if (mp_ID->GetFlipOutX())
1699  {
1700  // Loop over Z
1701  for (INTNB z=0; z<mp_ID->GetNbVoxZ(); z++)
1702  {
1703  // For efficiency
1704  INTNB base_z = z * mp_ID->GetNbVoxXY();
1705  // Loop over Y
1706  for (INTNB y=0; y<mp_ID->GetNbVoxY(); y++)
1707  {
1708  // For efficiency
1709  INTNB base_y = y * mp_ID->GetNbVoxX();
1710  // Half loop over X
1711  for (INTNB x_1=0; x_1<mp_ID->GetNbVoxX()/2; x_1++)
1712  {
1713  // Compute opposite X
1714  INTNB x_2 = mp_ID->GetNbVoxX() - 1 - x_1;
1715  // Compute both indices
1716  INTNB indice_1 = base_z + base_y + x_1;
1717  INTNB indice_2 = base_z + base_y + x_2;
1718  // Switch voxels
1719 /*
1720  FLTNB buffer = m4p_outputImage[fr][rg][cg][indice_1];
1721  m4p_outputImage[fr][rg][cg][indice_1] = m4p_outputImage[fr][rg][cg][indice_2];
1722  m4p_outputImage[fr][rg][cg][indice_2] = buffer;
1723 */
1724  FLTNB buffer = m4p_forwardImage[tbf][rbf][cbf][indice_1];
1725  m4p_forwardImage[tbf][rbf][cbf][indice_1] = m4p_forwardImage[tbf][rbf][cbf][indice_2];
1726  m4p_forwardImage[tbf][rbf][cbf][indice_2] = buffer;
1727  }
1728  }
1729  }
1730  }
1731  }
1732  }
1733  }
1734 
1735  // Normal end
1736  return 0;
1737 }
1738 
1739 // =====================================================================
1740 // ---------------------------------------------------------------------
1741 // ---------------------------------------------------------------------
1742 // =====================================================================
1743 
1745 {
1746  // Note: The output image is copied in the forward image at this step
1747 
1748  // 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
1749  if ( mp_ID->GetFOVOutPercent()<=0. &&
1750  mp_ID->GetNbSliceOutMask()==0 ) return 0;
1751  // Verbose
1752  if (m_verbose>=2)
1753  {
1754  Cout("oImageSpace::ApplyOutputFOVMasking() -> Mask output image" << endl);
1755  if (mp_ID->GetFOVOutPercent()>0.) Cout(" --> Mask transaxial FOV outside " << mp_ID->GetFOVOutPercent() << " %" << endl);
1756  if (mp_ID->GetNbSliceOutMask()>0) Cout(" --> Mask " << mp_ID->GetNbSliceOutMask() << " slices from both axial FOV limits" << endl);
1757  }
1758  // -----------------------------------------------
1759  // Transaxial FOV masking
1760  // -----------------------------------------------
1761  if (mp_ID->GetFOVOutPercent()>0.)
1762  {
1763  // Precast half the number of voxels over X and Y minus 1 (for efficiency)
1764  FLTNB flt_base_x = 0.5*((FLTNB)(mp_ID->GetNbVoxX()-1));
1765  FLTNB flt_base_y = 0.5*((FLTNB)(mp_ID->GetNbVoxY()-1));
1766  // Compute FOV elipse radius over X and Y, then squared
1767  FLTNB squared_radius_x = 0.5 * ((FLTNB)(mp_ID->GetNbVoxX())) * mp_ID->GetVoxSizeX()
1768  * mp_ID->GetFOVOutPercent() / 100.;
1769  squared_radius_x *= squared_radius_x;
1770  FLTNB squared_radius_y = 0.5 * ((FLTNB)(mp_ID->GetNbVoxY())) * mp_ID->GetVoxSizeY()
1771  * mp_ID->GetFOVOutPercent() / 100.;
1772  squared_radius_y *= squared_radius_y;
1773  // We assume that the computation of the distance from the center for a given
1774  // voxel and comparing it with the output FOV percentage costs more than performing
1775  // the loops in an inverse order compared to how the image is stored in memory.
1776  // Thus we begin the loops over X, then Y, then we test and if test passes, we
1777  // do the remaining loop over Z and over all dynamic dimensions.
1778  int x;
1779  #pragma omp parallel for private(x) schedule(guided)
1780  for (x=0; x<mp_ID->GetNbVoxX(); x++)
1781  {
1782  // Compute X distance from image center, then squared
1783  FLTNB squared_distance_x = (((FLTNB)x)-flt_base_x) * mp_ID->GetVoxSizeX();
1784  squared_distance_x *= squared_distance_x;
1785  // Loop over Y
1786  for (int y=0; y<mp_ID->GetNbVoxY(); y++)
1787  {
1788  // Compute Y distance from image center, then squared
1789  FLTNB squared_distance_y = (((FLTNB)y)-flt_base_y) * mp_ID->GetVoxSizeY();
1790  squared_distance_y *= squared_distance_y;
1791  // Test if the voxel is inside the FOV elipse, then we skip this voxel
1792  if ( squared_distance_x/squared_radius_x + squared_distance_y/squared_radius_y <= 1. ) continue;
1793  // Loop over Z
1794  for (int z=0; z<mp_ID->GetNbVoxZ(); z++)
1795  {
1796  // Compute global voxel index
1797  INTNB index = z*mp_ID->GetNbVoxXY() + y*mp_ID->GetNbVoxX() + x;
1798 /*
1799  // First loop on frames
1800  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1801  // Second loop on respiratory gates
1802  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1803  // Third loop on cardiac gates
1804  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1805 */
1806  // First loop on time basis functions
1807  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1808  // Second loop on respiratory basis functions
1809  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1810  // Third loop on cardiac basis functions
1811  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1812  // Put 0 into this voxel
1813 // m4p_outputImage[fr][rg][cg][index] = 0.;
1814  m4p_forwardImage[tbf][rbf][cbf][index] = 0.;
1815  }
1816  }
1817  }
1818  }
1819  // -----------------------------------------------
1820  // Axial FOV masking
1821  // -----------------------------------------------
1822  if (mp_ID->GetNbSliceOutMask()>0)
1823  {
1824  INTNB removed_slices = mp_ID->GetNbSliceOutMask();
1825 /*
1826  // First loop on frames
1827  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1828  {
1829  // Second loop on respiratory gates
1830  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1831  {
1832  // Third loop on cardiac gates
1833  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1834  {
1835 */
1836  // First loop on time basis functions
1837  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1838  {
1839  // Second loop on respiratory basis functions
1840  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1841  {
1842  // Third loop on cardiac basis functions
1843  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1844  {
1845  // Mask slices
1846  for (int z=0; z<removed_slices; z++)
1847  {
1848  // First slices
1849  INTNB base_z_first = z*mp_ID->GetNbVoxXY();
1850  // Loop over Y and X
1851  for (int i=0; i<mp_ID->GetNbVoxXY(); i++)
1852  {
1853  INTNB index = base_z_first + i;
1854 // m4p_outputImage[fr][rg][cg][index] = 0.;
1855  m4p_forwardImage[tbf][rbf][cbf][index] = 0.;
1856  }
1857  // Last slices
1858  INTNB base_z_last = (mp_ID->GetNbVoxZ()-1-z)*mp_ID->GetNbVoxXY();
1859  // Loop over Y and X
1860  for (int i=0; i<mp_ID->GetNbVoxXY(); i++)
1861  {
1862  INTNB index = base_z_last + i;
1863 // m4p_outputImage[fr][rg][cg][index] = 0.;
1864  m4p_forwardImage[tbf][rbf][cbf][index] = 0.;
1865  }
1866  }
1867  }
1868  }
1869  }
1870  }
1871  // End
1872  return 0;
1873 }
1874 
1875 // =====================================================================
1876 // ---------------------------------------------------------------------
1877 // ---------------------------------------------------------------------
1878 // =====================================================================
1879 
1881 {
1882  // If no mask image return
1883  if ( !IsLoadedMask() ) return 0;
1884  // Verbose
1885  if (m_verbose>=2)
1886  {
1887  Cout("oImageSpace::ApplyOutputMaskImage() -> Mask output image with the provided input mask image" << endl);
1888  }
1889  // -----------------------------------------------
1890  // Masking
1891  // -----------------------------------------------
1892  int v;
1893  #pragma omp parallel for private(v) schedule(guided)
1894  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1895  {
1896  // First loop on frames
1897  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1898  // Second loop on respiratory gates
1899  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1900  // Third loop on cardiac gates
1901  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1902  {
1903  if (!(mp_maskImage[v]>0.)) m4p_outputImage[fr][rg][cg][v] = 0.;
1904  }
1905  }
1906 
1907  // End
1908  return 0;
1909 }
1910 
1911 // =====================================================================
1912 // ---------------------------------------------------------------------
1913 // ---------------------------------------------------------------------
1914 // =====================================================================
1915 
1917 {
1918  // If no mask image return
1919  if ( !IsLoadedMask() ) return 0;
1920 
1921  // Verbose
1922  if (m_verbose>=2)
1923  {
1924  Cout("oImageSpace::ApplyMaskToSensitivity() -> Mask the sensitivity image with the provided input mask image" << endl);
1925  }
1926  // -----------------------------------------------
1927  // Masking
1928  // -----------------------------------------------
1929  int thread_0 = 0;
1930  int v;
1931  #pragma omp parallel for private(v) schedule(guided)
1932  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1933  {
1934  // First loop on frames
1935  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1936  // Second loop on respiratory gates
1937  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1938  // Third loop on cardiac gates
1939  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1940  {
1941  if (!(mp_maskImage[v]>0.)) m5p_sensitivity[thread_0][fr][rg][cg][v] = 0.;
1942  }
1943  }
1944 
1945  // End
1946  return 0;
1947 }
1948 
1949 // =====================================================================
1950 // ---------------------------------------------------------------------
1951 // ---------------------------------------------------------------------
1952 // =====================================================================
1953 
1954 int oImageSpace::ApplyMaskToBackwardImage(int a_imageIndex, int a_timeIndex, int a_respIndex, int a_cardIndex)
1955 {
1956  // If no mask image return
1957  if ( !IsLoadedMask() ) return 0;
1958 
1959  // Verbose
1960  if (m_verbose>=2)
1961  {
1962  Cout("oImageSpace::ApplyMaskToBackwardImage() -> Mask the backward image with the provided input mask image" << endl);
1963  }
1964  // -----------------------------------------------
1965  // Masking
1966  // -----------------------------------------------
1967  int thread_0 = 0;
1968  int v;
1969  #pragma omp parallel for private(v) schedule(guided)
1970  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1971  {
1972  // First loop on frames
1973  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1974  // Second loop on respiratory gates
1975  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1976  // Third loop on cardiac gates
1977  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1978  {
1979  if (!(mp_maskImage[v]>0.)) m6p_backwardImage[a_imageIndex][thread_0][a_timeIndex][a_respIndex][a_cardIndex][v] = 0.;
1980  }
1981  }
1982 
1983  // End
1984  return 0;
1985 }
1986 
1987 // =====================================================================
1988 // ---------------------------------------------------------------------
1989 // ---------------------------------------------------------------------
1990 // =====================================================================
1991 /*
1992  \fn SaveOutputImage
1993  \param a_iteration : current iteration index
1994  \param a_subset : current number of subsets (or -1 by default)
1995  \brief Call the interfile function to write output image on disk
1996  \return 0 if success, positive value otherwise
1997 */
1998 int oImageSpace::SaveOutputImage(int a_iteration, int a_subset)
1999 {
2000  if (m_verbose>=3) Cout("oImageSpace::SaveOutputImage ..."<< endl);
2001 
2002  // Get the output manager
2003  sOutputManager* p_output_manager = sOutputManager::GetInstance();
2004 
2005  // Get base name
2006  string path_to_img = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
2007 
2008  // Add a suffix for iteration
2009  if (a_iteration >= 0)
2010  {
2011  stringstream ss; ss << a_iteration + 1;
2012  path_to_img.append("_it").append(ss.str());
2013  }
2014 
2015  // Add a suffix for subset (if not negative by default), this means that we save a 'subset' image
2016  if (a_subset >= 0)
2017  {
2018  stringstream ss; ss << a_subset + 1;
2019  path_to_img.append("_ss").append(ss.str());
2020  }
2021 
2022  // We need one "mode" parameter which indicates if we would like to write 1 image file or several
2023  // Adjust functions according to that
2024  if (IntfWriteImgFile(path_to_img, m4p_outputImage, mp_ID, m_verbose) )
2025  {
2026  Cerr("***** oImageSpace::SaveOutputImage()-> Error writing Interfile of output image !" << endl);
2027  return 1;
2028  }
2029 
2030  // Normal end
2031  return 0;
2032 }
2033 
2034 // =====================================================================
2035 // ---------------------------------------------------------------------
2036 // ---------------------------------------------------------------------
2037 // =====================================================================
2038 
2039 int oImageSpace::SaveOutputBasisCoefficientImage(int a_iteration, int a_subset)
2040 {
2041  if (m_verbose>=3) Cout("oImageSpace::SaveOutputBasisCoefficientImage ..."<< endl);
2042 
2043  // Get the output manager
2044  sOutputManager* p_output_manager = sOutputManager::GetInstance();
2045 
2046  // Get base name
2047  string path_to_img = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
2048 
2049  // Add a suffix for iteration
2050  if (a_iteration >= 0)
2051  {
2052  stringstream ss; ss << a_iteration + 1;
2053  path_to_img.append("_it").append(ss.str());
2054  }
2055 
2056  // Add a suffix for subset (if not negative by default), this means that we save a 'subset' image
2057  if (a_subset >= 0)
2058  {
2059  stringstream ss; ss << a_subset + 1;
2060  path_to_img.append("_ss").append(ss.str());
2061  }
2062 
2063  // At this step, the m4p_forwardImage contains the actual image basis coefficients for output
2065  {
2066  Cerr("***** oImageSpace::SaveOutputImage()-> Error writing Interfile of output image !" << endl);
2067  return 1;
2068  }
2069 
2070  // Normal end
2071  return 0;
2072 }
2073 
2074 
2075 
2076 // =====================================================================
2077 // ---------------------------------------------------------------------
2078 // ---------------------------------------------------------------------
2079 // =====================================================================
2080 /*
2081  \fn SaveDebugImage
2082  \param a_name : output name of the image
2083  \brief Just a debug function dedicated to write any kind of image on disk in raw format, for debugging purposes
2084 */
2085 void oImageSpace::SaveDebugImage(const string& a_name)
2086 {
2087  #ifdef CASTOR_DEBUG
2088  if(m_verbose>=3) Cout("oImageSpace::SaveDebugImage ..."<< endl);
2089  #endif
2090 
2091  ofstream output_file;
2092  output_file.open(a_name.c_str(), ios::binary | ios::out);
2093 
2094  //for (int fr=0; fr<ap_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
2095  // for (int rg=0; rg<ap_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
2096  // for (int cg=0; cg<ap_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
2097  // {
2098  // Write file
2099  output_file.write(reinterpret_cast<char*>(m6p_backwardImage[0][0][0][0][0]), mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
2100  // Close file
2101  output_file.close();
2102  // }
2103  // }
2104  //}
2105 }
2106 
2107 
2108 
2109 
2110 // =====================================================================
2111 // ---------------------------------------------------------------------
2112 // ---------------------------------------------------------------------
2113 // =====================================================================
2114 /*
2115  \fn PrepareForwardImage
2116  \brief Copy current image matrix in the forward-image buffer matrix
2117 */
2119 {
2120  if(m_verbose>=3) Cout("oImageSpace::PrepareForwardImage ..."<< endl);
2121 
2122  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2123  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2124  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2125  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2126  m4p_forwardImage[tbf][rbf][cbf][v] = m4p_image[tbf][rbf][cbf][v];
2127 }
2128 
2129 
2130 
2131 
2132 // =====================================================================
2133 // ---------------------------------------------------------------------
2134 // ---------------------------------------------------------------------
2135 // =====================================================================
2136 /*
2137  \fn Reduce
2138  \brief Merge parallel results into the matrix of the backward image matrix of the first thread. Also for MPI.
2139  \todo Choose computation alternative (do not know yet which one is faster...)
2140 */
2142 {
2143  #if defined(CASTOR_OMP) || defined(CASTOR_MPI)
2144  if (m_verbose>=3) Cout("oImageSpace::Reduce() -> Merge parallel results" << endl);
2145  #endif
2146 
2147  // --------------------------------------------------------------------------------
2148  // Step 1: merge multi-threads results
2149  // --------------------------------------------------------------------------------
2150 
2151  #ifdef CASTOR_OMP
2152  if (m_verbose>=3) Cout(" --> Over threads ..." << endl);
2153 
2154  // Special case here where it appears to be always beneficial to use as many threads
2155  // as the number of threads used for projections. So we set it here, and set it back
2156  // at the end.
2157  omp_set_num_threads(mp_ID->GetNbThreadsForProjection());
2158 
2159  // todo : Choose computation alternative (do not know yet which one is faster...)
2160  int alternative = 2;
2161 
2162  // Alternative 1: standard loops based on the pointers hierarchy but mono-thread
2163  if (alternative==1 || mp_ID->GetNbThreadsForImageComputation()==1)
2164  {
2165  for (int img=0; img<m_nbBackwardImages; img++)
2166  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
2167  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2168  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2169  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2170  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2171  m6p_backwardImage[img][0][tbf][rbf][cbf][v] += m6p_backwardImage[img][th][tbf][rbf][cbf][v];
2172  }
2173  // Alternative 2: multi-thread merging with the voxels loop at first to be thread safe
2174  // Maybe faster than alternative 1, but the first loop on voxels breaks the memory order... have to try
2175  else if (alternative==2)
2176  {
2177  int v;
2178  #pragma omp parallel for private(v) schedule(guided)
2179  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2180  {
2181  for (int img=0; img<m_nbBackwardImages; img++)
2182  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
2183  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2184  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2185  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2186  m6p_backwardImage[img][0][tbf][rbf][cbf][v] += m6p_backwardImage[img][th][tbf][rbf][cbf][v];
2187  }
2188  }
2189  // If on-the-fly sensitivity, then do it also for it.
2190  if (!m_loadedSensitivity)
2191  {
2192  if (alternative==1 || mp_ID->GetNbThreadsForImageComputation()==1)
2193  {
2194  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
2195  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2196  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
2197  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
2198  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2199  m5p_sensitivity[0][fr][rg][cg][v] += m5p_sensitivity[th][fr][rg][cg][v];
2200  }
2201  else if (alternative==2)
2202  {
2203  int v;
2204  #pragma omp parallel for private(v) schedule(guided)
2205  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2206  {
2207  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
2208  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2209  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
2210  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
2211  m5p_sensitivity[0][fr][rg][cg][v] += m5p_sensitivity[th][fr][rg][cg][v];
2212  }
2213  }
2214  }
2215  // Set back the number of threads to the one for image computation
2216  omp_set_num_threads(mp_ID->GetNbThreadsForImageComputation());
2217  #endif
2218 
2219  // --------------------------------------------------------------------------------
2220  // Step 2: merge multi-instance results (MPI)
2221  // --------------------------------------------------------------------------------
2222 
2223  #ifdef CASTOR_MPI
2224  // We use the MPI_IN_PLACE as the send buffer so that it is the same as the result buffer
2225  if (m_verbose>=3) Cout(" --> Over instances ..." << endl);
2226 
2227  // Merge backward images
2228  for (int img=0; img<m_nbBackwardImages; img++)
2229  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2230  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2231  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2232  MPI_Allreduce( MPI_IN_PLACE, &m6p_backwardImage[img][0][tbf][rbf][cbf][0], mp_ID->GetNbVoxXYZ(), FLTNBMPI, MPI_SUM, MPI_COMM_WORLD );
2233 
2234  // If on-the-fly sensitivity, then do it also for it
2235  if (!m_loadedSensitivity)
2236  {
2237  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2238  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
2239  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
2240  MPI_Allreduce( MPI_IN_PLACE, &m5p_sensitivity[0][fr][rg][cg][0], mp_ID->GetNbVoxXYZ(), FLTNBMPI, MPI_SUM, MPI_COMM_WORLD );
2241  }
2242  #endif
2243 }
2244 
2245 
2246 
2247 
2248 
2249 // =====================================================================
2250 // ---------------------------------------------------------------------
2251 // ---------------------------------------------------------------------
2252 // =====================================================================
2253 /*
2254  \fn CleanNeverVisitedVoxels
2255  \brief Based on the visitedVoxelsImage, clean the never visited voxels in the image.
2256  This function must be called at the end of each iteration
2257 */
2259 {
2260  if(m_verbose>=3) Cout("oImageSpace::CleanNeverVisitedVoxels ..."<< endl);
2261  // Multi-threaded loop over voxels
2262  INTNB v;
2263  #pragma omp parallel for private(v) schedule(guided)
2264  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2265  {
2266  // Clean this image voxel if the voxel was never visited
2267  if (mp_visitedVoxelsImage[v]==0.)
2268  {
2269  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2270  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2271  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2272  m4p_image[tbf][rbf][cbf][v] = 0.;
2273  }
2274  // Reset the image
2275  mp_visitedVoxelsImage[v] = 0.;
2276  }
2277 }
2278 
2279 
2280 
2281 
2282 
2283 
2284 
2285 
2286 
2287 // LIST-MODE SENSITIVITY GENERATION FUNCTIONS
2288 // =====================================================================
2289 // ---------------------------------------------------------------------
2290 // ---------------------------------------------------------------------
2291 // =====================================================================
2292 /*
2293  \fn LMS_InstantiateImage()
2294  \brief Allocate memory for the main image matrices (for list-mode sensitivity generation)
2295  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2296 */
2298 {
2299  if(m_verbose>=3) Cout("oImageSpace::LMS_InstantiateImage ..."<< endl);
2300 
2301  m4p_image = new FLTNB***[mp_ID->GetNbTimeFrames()];
2302 
2303  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2304  {
2305  m4p_image[fr] = new FLTNB**[mp_ID->GetNb1stMotImgsForLMS(fr)];
2306 
2307  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2308  {
2309  m4p_image[fr][rg] = new FLTNB*[mp_ID->GetNb2ndMotImgsForLMS()];
2310 
2311  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2312  m4p_image[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
2313  }
2314  }
2315 }
2316 
2317 
2318 
2319 
2320 // =====================================================================
2321 // ---------------------------------------------------------------------
2322 // ---------------------------------------------------------------------
2323 // =====================================================================
2324 /*
2325  \fn LMS_DeallocateImage()
2326  \brief Free memory for the main image matrices (for list-mode sensitivity generation)
2327  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2328 */
2330 {
2331  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateImage ..."<< endl);
2332 
2333  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2334  {
2335  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2336  {
2337  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2338  if(m4p_image[fr][rg][cg] != NULL) delete m4p_image[fr][rg][cg];
2339 
2340  if(m4p_image[fr][rg] != NULL) delete[] m4p_image[fr][rg];
2341  }
2342  if(m4p_image[fr] != NULL) delete[] m4p_image[fr] ;
2343  }
2344 
2345  if(m4p_image != NULL) delete[] m4p_image;
2346  m4p_image = NULL;
2347 }
2348 
2349 
2350 
2351 // =====================================================================
2352 // ---------------------------------------------------------------------
2353 // ---------------------------------------------------------------------
2354 // =====================================================================
2355 /*
2356  \fn LMS_InstantiateForwardImage()
2357  \brief Allocate memory for the forward image matrices (for list-mode sensitivity generation)
2358  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2359 */
2361 {
2362  if(m_verbose>=3) Cout("oImageSpace::LMS_InstantiateForwardImage ..."<< endl);
2363 
2365 
2366  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2367  {
2369 
2370  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2371  {
2373 
2374  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2375  {
2376  m4p_forwardImage[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
2377  }
2378  }
2379 
2380  }
2381 }
2382 
2383 
2384 
2385 // =====================================================================
2386 // ---------------------------------------------------------------------
2387 // ---------------------------------------------------------------------
2388 // =====================================================================
2389 /*
2390  \fn LMS_DeallocateForwardImage()
2391  \brief Free memory for the forward image matrices (for list-mode sensitivity generation)
2392  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2393 */
2395 {
2396  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateForwardImage ..."<< endl);
2397 
2398  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2399  {
2400  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2401  {
2402  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2403  {
2404  if(m4p_forwardImage[fr][rg][cg] != NULL) delete m4p_forwardImage[fr][rg][cg];
2405  }
2406 
2407  if(m4p_forwardImage[fr][rg] != NULL) delete[] m4p_forwardImage[fr][rg];
2408  }
2409 
2410  if(m4p_forwardImage[fr] != NULL) delete[] m4p_forwardImage[fr];
2411  }
2412 
2413  if(m4p_forwardImage != NULL) delete[] m4p_forwardImage;
2414  m4p_forwardImage = NULL;
2415 }
2416 
2417 // =====================================================================
2418 // ---------------------------------------------------------------------
2419 // ---------------------------------------------------------------------
2420 // =====================================================================
2421 /*
2422  \fn LMS_InstantiateSensitivityImage()
2423  \brief Allocate memory for the sensitivity image matrices (for list-mode sensitivity generation)
2424  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2425 */
2427 {
2428  if(m_verbose>=3) Cout("oImageSpace::LMS_InstantiateSensitivityImage ..."<< endl);
2429 
2430  m5p_sensitivity = new FLTNB****[1];
2431  m5p_sensitivity[0] = new FLTNB***[mp_ID->GetNbTimeFrames()];
2432  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2433  {
2434  m5p_sensitivity[0][fr] = new FLTNB**[mp_ID->GetNb1stMotImgsForLMS(fr)];
2435  for (int rg=0; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
2436  {
2437  m5p_sensitivity[0][fr][rg] = new FLTNB*[mp_ID->GetNb2ndMotImgsForLMS()];
2438  for (int cg=0; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
2439  m5p_sensitivity[0][fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
2440  }
2441  }
2442 }
2443 
2444 
2445 
2446 
2447 // =====================================================================
2448 // ---------------------------------------------------------------------
2449 // ---------------------------------------------------------------------
2450 // =====================================================================
2451 /*
2452  \fn LMS_DeallocateSensitivityImage()
2453  \brief Free memory for the sensitivity image matrices (for list-mode sensitivity generation)
2454  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2455 */
2457 {
2458  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateSensitivityImage ..."<< endl);
2459 
2460  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2461  {
2462  for (int rg=0; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
2463  {
2464  for (int cg=0; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
2465  {
2466  if(m5p_sensitivity[0][fr][rg][cg] != NULL) delete m5p_sensitivity[0][fr][rg][cg];
2467  }
2468  if(m5p_sensitivity[0][fr][rg] != NULL)delete m5p_sensitivity[0][fr][rg];
2469  }
2470  if(m5p_sensitivity[0][fr] != NULL)delete m5p_sensitivity[0][fr];
2471  }
2472  if(m5p_sensitivity[0] != NULL) delete[] m5p_sensitivity[0];
2473 
2474  if(m5p_sensitivity != NULL) delete[] m5p_sensitivity;
2475  m5p_sensitivity = NULL;
2476 }
2477 
2478 
2479 
2480 
2481 // =====================================================================
2482 // ---------------------------------------------------------------------
2483 // ---------------------------------------------------------------------
2484 // =====================================================================
2485 /*
2486  \fn LMS_CopyAtnToImage()
2487  \brief Copy the attenuation image contained in the 'm2p_attenuation' matrix inside the m2p_image matrix (for list-mode sensitivity generation)
2488  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2489 */
2491 {
2492  if(m_verbose>=3) Cout("oImageSpace::LMS_CopyAtnToImage ..."<< endl);
2493 
2494  if(m4p_attenuation != NULL)
2495  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2496  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2497  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2498  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2499  {
2500  m4p_image[fr][rg][cg][v] = m4p_attenuation[fr][rg][cg][v];
2501  }
2502  else
2503  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2504  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2505  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2506  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2507  {
2508  m4p_image[fr][rg][cg][v] = 0;
2509  }
2510 
2511 }
2512 
2513 
2514 // =====================================================================
2515 // ---------------------------------------------------------------------
2516 // ---------------------------------------------------------------------
2517 // =====================================================================
2518 /*
2519  \fn LMS_CopyAtnToForwardImage()
2520  \param a_use1stMotion
2521  \param a_use2ndMotion
2522  \brief Copy the attenuation image contained in the 'm2p_attenuation' matrix inside the m2p_forwardImage matrix (for list-mode sensitivity generation)
2523  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2524 */
2525 void oImageSpace::LMS_CopyAtnToForwardImage(bool a_use1stMotion, bool a_use2ndMotion)
2526 {
2527  if(m_verbose>=3) Cout("oImageSpace::LMS_CopyAtnToForwardImage ..."<< endl);
2528 
2529  if(m4p_attenuation != NULL)
2530  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2531  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2532  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2533  {
2534  // If motion is enabled, copy the attenuation image in all forward image matrices
2535  int rg_atn = a_use1stMotion ? 0 : rg ;
2536  int cg_atn = a_use1stMotion ? 0 : cg ;
2537 
2538  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2539  {
2540  m4p_forwardImage[fr][rg][cg][v] = m4p_attenuation[fr][rg_atn][cg_atn][v];
2541  }
2542  }
2543  else
2544  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2545  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2546  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2547  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2548  {
2549  m4p_forwardImage[fr][rg][cg][v] = 0;
2550  }
2551 }
2552 
2553 // =====================================================================
2554 // ---------------------------------------------------------------------
2555 // ---------------------------------------------------------------------
2556 // =====================================================================
2557 /*
2558  \fn LMS_CopyBackwardToSensitivityImage()
2559  \brief Copy the backward image containing the result of the sensitivity back-projection, into the sensitivity matrix.
2560  \details This function is dedicated to list-mode sensitivity (LMS) generation, it is useful to apply convolution and to save the image.
2561 */
2563 {
2564  if(m_verbose>=3) Cout("oImageSpace::LMS_CopyBackwardToSensitivity ..."<< endl);
2565  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2566  for (int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2567  for (int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2568  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2569  m5p_sensitivity[0][fr][rg][cg][v] = m6p_backwardImage[0][0][fr][rg][cg][v];
2570 }
2571 
2572 
2573 // =====================================================================
2574 // ---------------------------------------------------------------------
2575 // ---------------------------------------------------------------------
2576 // =====================================================================
2577 /*
2578  \fn LMS_PrepareForwardImage()
2579  \brief Copy current image in forward-image buffer (for list-mode sensitivity generation)
2580  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2581 */
2583 {
2584  if(m_verbose>=3) Cout("oImageSpace::LMS_PrepareForwardImage ..."<< endl);
2585 
2586  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2587  for (int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2588  for (int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2589  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2590  m4p_forwardImage[fr][rg][cg][v] = m4p_image[fr][rg][cg][v];
2591 }
2592 
2593 
2594 
2595 
2596 // =====================================================================
2597 // ---------------------------------------------------------------------
2598 // ---------------------------------------------------------------------
2599 // =====================================================================
2600 
2601 void oImageSpace::ReduceBackwardImage(int a_imageIndex, int a_timeIndex, int a_respIndex, int a_cardIndex)
2602 {
2603  #if defined(CASTOR_OMP) || defined(CASTOR_MPI)
2604  if (m_verbose>=3) Cout("oImageSpace::ReduceBackwardImage() -> Merge parallel results" << endl);
2605  #endif
2606 
2607  // --------------------------------------------------------------------------------
2608  // Step 1: merge multi-threads results
2609  // --------------------------------------------------------------------------------
2610 
2611  #ifdef CASTOR_OMP
2612  // Verbose
2613  if (m_verbose>=3) Cout(" --> Over threads ..." << endl);
2614  // Do it
2615  for (int th=1 ; th<mp_ID->GetNbThreadsForProjection() ; th++)
2616  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2617  m6p_backwardImage[a_imageIndex][0][a_timeIndex][a_respIndex][a_cardIndex][v] += m6p_backwardImage[a_imageIndex][th][a_timeIndex][a_respIndex][a_cardIndex][v];
2618  #endif
2619 
2620  // --------------------------------------------------------------------------------
2621  // Step 2: merge multi-instance results (MPI)
2622  // --------------------------------------------------------------------------------
2623 
2624  #ifdef CASTOR_MPI
2625  // We use the MPI_IN_PLACE as the send buffer so that it is the same as the result buffer
2626  if (m_verbose>=3) Cout(" --> Over instances ..." << endl);
2627  // Merge backward images
2628  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 );
2629  #endif
2630 }
2631 
2632 
2633 
2634 
2635 // =====================================================================
2636 // ---------------------------------------------------------------------
2637 // ---------------------------------------------------------------------
2638 // =====================================================================
2639 /*
2640  \fn LMS_SaveSensitivityImage
2641  \param a_pathToSensitivityImage : path to the sensitivity image (should be provided only in list-mode reconstruction)
2642  \param ap_DeformationManager : Pointer to the deformation manager objet (required to retrieve the number of gates in the sensitivity image)
2643  \brief Call the interfile function to write the sensitivity image on disk
2644  \details If image deformation is enabled for respiratory/cardiac gated data, the gated images are summed up into one image and normalize
2645  \return 0 if success, positive value otherwise
2646  \todo Interfile management : if the number of sensitivity image to write is different to the actual number of image
2647  as it could be the case for dynamic imaging, if one sensitivity image is required for the whole dynamic serie,
2648  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)
2649 */
2650 int oImageSpace::LMS_SaveSensitivityImage(const string& a_pathToSensitivityImage, oDeformationManager* ap_DeformationManager)
2651 {
2652  if(m_verbose>=3) Cout("oImageSpace::LMS_SaveSensitivityImage ..."<< endl);
2653 
2654  // If image-based motion is enabled, list-mode sensitivity images code generates a set of several provisional images of the reference position.
2655  // Each image is generated using one different forward/backward deformations parameters which will be used during reconstruction.
2656  // - For gated (respiratory/cardiac) motion correction, a simple average of these images is performed here for each frame.
2657  // - For timestamp-based motion (involuntary patient motion), the averaging depends on which transformation occurs in each frame.
2658  // The provisional sensitivity images are weighs accordingly to the duration of each motion subset inside the frame.
2659  // (i.e the sensitivity image of a 10 min frame, which contains a transformation at 9mn, will be composed of an average of two sensitivity images with weights of 0.9 and 0.1)
2660 
2661  int nb_reco_card_images = ap_DeformationManager->GetNbSensImagesCardDeformation(mp_ID->GetNb2ndMotImgsForLMS());
2662 
2663  // Average of the cardiac images, if cardiac gating is enabled (no entry in the loop of cardiac images otherwise)
2664  if (ap_DeformationManager->GetNbSensImagesCardDeformation(mp_ID->GetNb2ndMotImgsForLMS()) == 0)
2665  {
2666  nb_reco_card_images = 1;
2667  FLTNB card_normalization = mp_ID->GetNb2ndMotImgsForLMS();
2668 
2669  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2670  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames(); fr++)
2671  for (int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2672  {
2673  for (int cg=1 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2674  m5p_sensitivity[0][fr][rg][0][v] += m5p_sensitivity[0][fr][rg][cg][v];
2675  m5p_sensitivity[0][fr][rg][0][v] /= card_normalization;
2676  }
2677  }
2678 
2679  // Average of the respiratory images, if respiratory gating is enabled (no entry in the loop of respiratory images)
2680  if (!mp_ID->IsPMotionEnabled() // No patient motion correction
2681  && ap_DeformationManager->GetNbSensImagesRespDeformation(mp_ID->GetNb1stMotImgsForLMS(0)) == 0)
2682  {
2683  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames(); fr++)
2684  {
2685  FLTNB resp_normalization = mp_ID->GetNb1stMotImgsForLMS(fr);
2686 
2687  for (int cg=0 ; cg<nb_reco_card_images ; cg++) // Dual gating : Be sure to get the right number of card images
2688  {
2689  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2690  {
2691  for (int rg=1 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2692  m5p_sensitivity[0][fr][0][cg][v] += m5p_sensitivity[0][fr][rg][cg][v];
2693 
2694  m5p_sensitivity[0][fr][0][cg][v] /= resp_normalization;
2695  }
2696  }
2697  }
2698  }
2699 
2700  // Average of the sensitivity images when involuntary patient motion correction is enabled
2701  if (mp_ID->IsPMotionEnabled() )
2702  {
2703  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames(); fr++)
2704  {
2705  //FLTNB ipm_normalization = GetNb1stMotImgsForLMS(fr);
2706 
2707  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2708  {
2709 
2710 
2711  //for (int ipm=1 ; ipm<mp_ID->GetNb1stMotImgsForLMS(fr) ; ipm++)
2712  // m5p_sensitivity[0][fr][0][0][v] += m5p_sensitivity[0][fr][ipm][0][v] ;
2713 
2714  //m5p_sensitivity[0][fr][0][0][v] /= ipm_normalization;
2715 
2716  FLTNB sensitivity_value_avg = 0.;
2717 
2718  // Add the contribution of the sensitivity images with appropriate weighing
2719  for (int ipm=0 ; ipm<mp_ID->GetNb1stMotImgsForLMS(fr) ; ipm++)
2720  sensitivity_value_avg += m5p_sensitivity[0][fr][ipm][0][v] * mp_ID->GetListPMotionWeightInFrameForLMS(fr, ipm);
2721 
2722  // Recover the average value in the first sensitivity image
2723  m5p_sensitivity[0][fr][0][0][v] = sensitivity_value_avg;
2724  }
2725 
2726  }
2727  }
2728 
2729 
2730 
2731  if (SaveSensitivityImage(a_pathToSensitivityImage) )
2732  {
2733  Cerr("***** oImageSpace::LMS_SaveSensitivityImage()-> Error writing Sensitivity image !" << endl);
2734  return 1;
2735  }
2736 
2737  return 0;
2738 }
2739 
2740 
2741 
2742 
2743 // =====================================================================
2744 // ---------------------------------------------------------------------
2745 // ---------------------------------------------------------------------
2746 // =====================================================================
2747 /*
2748  \fn SaveSensitivityImage
2749  \param a_pathToSensitivityImage : path to the sensitivity image (should be provided only in list-mode reconstruction)
2750  \brief Call the interfile function to write the sensitivity image on disk
2751  \return 0 if success, positive value otherwise
2752 */
2753 int oImageSpace::SaveSensitivityImage(const string& a_pathToSensitivityImage)
2754 {
2755  if(m_verbose>=3) Cout("oImageSpace::SaveSensitivityImage ..."<< endl);
2756 
2757  if (IntfWriteImgFile(a_pathToSensitivityImage, m5p_sensitivity[0], mp_ID, m_verbose) )
2758  {
2759  Cerr("***** oImageSpace::SaveSensitivityImage()-> Error writing Interfile of output image !" << endl);
2760  return 1;
2761  }
2762 
2763  return 0;
2764 }
2765 
2766 
2767 
2768 
2769 // PROJECTION SCRIPT FUNCTIONS
2770 // =====================================================================
2771 // ---------------------------------------------------------------------
2772 // ---------------------------------------------------------------------
2773 // =====================================================================
2774 /*
2775  \fn PROJ_InstantiateProjectionImage()
2776  \param a_nbProjs : a number of projection slices in the projection
2777  \param a_nbPixels : a total number of pixels in the projection slices
2778  \brief Instanciate and initialize projection image for analytical projection
2779  \details This function is currently only dedicated to SPECT projection, and used by the analytical projection script
2780  \todo support for PET (requires projection format)
2781 */
2782 void oImageSpace::PROJ_InstantiateProjectionImage(int a_nbProjs, int a_nbPixels)
2783 {
2784  if(m_verbose>=3) Cout("oImageSpace::PROJ_InstantiateProjectionImage ..."<< endl);
2785 
2786  m2p_projectionImage = new FLTNB*[a_nbProjs];
2787 
2788  for(int p=0 ; p<a_nbProjs ; p++)
2789  {
2790  m2p_projectionImage[p] = new FLTNB[a_nbPixels];
2791  for(int px=0 ; px<a_nbPixels ; px++)
2792  m2p_projectionImage[p][px] = 0;
2793  }
2794 }
2795 
2796 
2797 
2798 
2799 // =====================================================================
2800 // ---------------------------------------------------------------------
2801 // ---------------------------------------------------------------------
2802 // =====================================================================
2803 /*
2804  \fn PROJ_DeallocateProjectionImage()
2805  \param a_nbProjs : a number of projection slices in the projection
2806  \brief Free memory for the projection image for analytical projection
2807  \details This function is currently only dedicated to SPECT projection, and used by the analytical projection script
2808  \todo support for PET (requires projection format)
2809 */
2811 {
2812  if(m_verbose>=3) Cout("oImageSpace::PROJ_DeallocateProjectionImage ..."<< endl);
2813 
2814  for(int p=0 ; p<a_nbProjs ; p++)
2815  {
2816  delete[] m2p_projectionImage[p];
2817  }
2818 
2819  delete[] m2p_projectionImage;
2820  m2p_projectionImage = NULL;
2821 }
2822 
2823 
2824 
2825 
2826 // =====================================================================
2827 // ---------------------------------------------------------------------
2828 // ---------------------------------------------------------------------
2829 // =====================================================================
2830 /*
2831  \fn PROJ_InitImage()
2832  \param a_pathToInitialImage : path to the image to project
2833  \brief Load the initial image for the analytical projection
2834  \return 0 if success, positive value otherwise
2835 */
2836 int oImageSpace::PROJ_InitImage(const string& a_pathToInitialImage)
2837 {
2838  if(m_verbose>=3) Cout("oImageSpace::PROJ_InitImage ..."<< endl);
2839 
2840  if (!a_pathToInitialImage.empty()) // Image initiale
2841  {
2842  if(PROJ_LoadInitialImage(a_pathToInitialImage) )
2843  {
2844  Cerr("***** oImageSpace::PROJ_InitImage()-> Error while trying to read file at " << a_pathToInitialImage << endl);
2845  return 1;
2846  }
2847  }
2848  else
2849  {
2850  {
2851  Cerr("***** oImageSpace::PROJ_InitImage()-> No projected image provided ! " << endl);
2852  return 1;
2853  }
2854  }
2855  return 0;
2856 }
2857 
2858 
2859 
2860 
2861 // =====================================================================
2862 // ---------------------------------------------------------------------
2863 // ---------------------------------------------------------------------
2864 // =====================================================================
2865 /*
2866  \fn PROJ_LoadInitialImage()
2867  \param a_pathToImage : path to the image to project
2868  \brief Load the initial image for the analytical projection
2869  \return 0 if success, positive value otherwise
2870 */
2871 int oImageSpace::PROJ_LoadInitialImage(const string& a_pathToImage)
2872 {
2873  if(m_verbose>=3) Cout("oImageSpace::PROJ_LoadInitialImage ..."<< endl);
2874 
2875  ifstream image_file;
2876  image_file.open(a_pathToImage.c_str(), ios::in | ios::binary);
2877 
2878  if(!image_file.is_open())
2879  {
2880  Cerr("***** oImageSpace::PROJ_LoadInitialImage()-> Error reading file !" << endl);
2881  return 1;
2882  }
2883  else
2884  {
2885  // Interfile image reading (INTF_LERP_DISABLED = no interpolation allowed)
2887  {
2888  Cerr("***** oImageSpace::PROJ_LoadInitialImage()-> Error reading Interfile : " << a_pathToImage << " !" << endl);
2889  return 1;
2890  }
2891 
2892  }
2893  image_file.close();
2894  return 0;
2895 }
2896 
2897 
2898 
2899 
2900 
2901 // =====================================================================
2902 // ---------------------------------------------------------------------
2903 // ---------------------------------------------------------------------
2904 // =====================================================================
2905 /*
2906  \fn PROJ_SaveProjectionImage
2907  \brief Save an image of the projected data
2908  \return 0 if success, positive value otherwise
2909  \todo Support for PET projeted data
2910  implement interfile file recovery in a scanner class function
2911 */
2913 {
2914  if(m_verbose>=3) Cout("oImageSpace::PROJ_SaveProjectionImage ..."<< endl);
2915 
2916  string img_name = "_ProjectionImage";
2917 
2918  // Initialize interfile fields structure
2919  Intf_fields Img_fields;
2920  IntfKeySetFieldsOutput(&Img_fields , mp_ID);
2921 
2922  Img_fields.process_status = "acquired";
2923 
2924  // Get specific info about projection
2926 
2927  // Recover the values in the interfile structure.
2928  // todo : This step for PET projection
2929  // todo : Implement this step in scanner class functions
2930  // todo : For spect, deal with systems with several detector heads (currently assume one head)
2931  if(p_ScanMgr->GetScannerType() == 2)
2932  {
2933  uint16_t nb_projs, nb_heads;
2934  FLTNB zoom;
2935  uint16_t nb_bins[2];
2936  FLTNB pix_size[2];
2937  FLTNB* p_angles;
2938  FLTNB* p_radius;
2939  int dir_rot;
2940  p_ScanMgr->GetSPECTSpecificParameters(&nb_projs,
2941  &nb_heads,
2942  &zoom,
2943  nb_bins,
2944  pix_size,
2945  p_angles,
2946  p_radius,
2947  &dir_rot);
2948 
2949  Img_fields.nb_projections = nb_projs;
2950  Img_fields.nb_detector_heads = nb_heads;
2951  Img_fields.mtx_size[0] = nb_bins[0];
2952  Img_fields.mtx_size[1] = nb_bins[1];
2953  Img_fields.vox_size[0] = pix_size[0];
2954  Img_fields.vox_size[1] = pix_size[1];
2955  Img_fields.extent_rotation = 360; // 360 by default
2956  // Check rotation direction from the two first angles
2957  Img_fields.direction_rotation = (dir_rot == GEO_ROT_CCW)?
2958  "CCW" :
2959  "CW";
2960  Img_fields.first_angle = p_angles[0];
2961 
2962  // Note: In Interfile v3.3, doesn't seem to have a field to provide
2963  // individually each angle and Center of rotation.
2964  // Looks like it was planned for ulterior version, at least for the
2965  // center of rotation
2966  // For now, we just write each projection angle and radius as an array key
2967 
2968  string angles_str = "{";
2969  string radius_str = "{";
2970 
2971  // Flag to check if the radius is similar for each projection angle
2972  bool has_single_radius = true;
2973 
2974  for(uint16_t p=0 ; p<nb_projs ; p++)
2975  {
2976  stringstream ss_a, ss_r;
2977  // projection angles
2978  ss_a << p_angles[p];
2979  angles_str.append(ss_a.str());
2980  (p<nb_projs-1) ? angles_str.append(",") : angles_str.append("}");
2981 
2982  // radius
2983  ss_r << p_radius[p];
2984  radius_str.append(ss_r.str());
2985  (p<nb_projs-1) ? radius_str.append(",") : radius_str.append("}");
2986  if(p_radius[p] != p_radius[0])
2987  has_single_radius = false;
2988 
2989  // Some interfile editors struggle with long line, perhaps because
2990  // they are stored in char[256] or something like that
2991  // Hence, break line after a certain number of elements
2992  if((p+1)%30 == 0)
2993  {
2994  angles_str.append("\n");
2995  radius_str.append("\n");
2996  }
2997  }
2998 
2999  // If there is one common radius, write its value in the related string
3000  if(has_single_radius)
3001  {
3002  stringstream ss;
3003  ss << p_radius[0];
3004  radius_str = ss.str();
3005  }
3006 
3007  Img_fields.projection_angles = angles_str;
3008  Img_fields.radius = radius_str;
3009 
3010 
3011  // Common code to all modalities
3012  sOutputManager* p_output_manager = sOutputManager::GetInstance();
3013  string img_file_name = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
3014  img_file_name.append("_ProjImage");
3015 
3016  if(IntfWriteProjFile(img_file_name, m2p_projectionImage, mp_ID, Img_fields, m_verbose) )
3017  {
3018  Cerr("***** oImageSpace::PROJ_SaveProjectionImage()-> Error writing Interfile of output image !" << endl);
3019  return 1;
3020  }
3021  }
3022  return 0;
3023 }
3024 
3025 
int ApplyMaskToBackwardImage(int a_imageIndex, int a_timeIndex, int a_respIndex, int a_cardIndex)
#define INTF_LERP_DISABLED
void IntfKeySetFieldsOutput(Intf_fields *ap_IF, oImageDimensionsAndQuantification *ap_ID)
Init the keys of the Interfile header of an image to be written on disk.
void LMS_InstantiateImage()
Allocate memory for the main image matrices (for list-mode sensitivity generation) ...
static sScannerManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
void DeallocateMiscellaneousImage()
Deallocate all allocated miscellaneous images.
void DeallocateBackwardImageFromDynamicBasis()
Free memory for the backward image matrices.
void LMS_CopyAtnToForwardImage(bool a_use1stMotion, bool a_use2ndMotion)
int GetNbCardBasisFunctions()
Get the number of cardiac basis functions.
#define INTF_LERP_ENABLED
int InitMaskImage(const string &a_pathToImage)
int PROJ_InitImage(const string &a_pathToInitialImage)
FLTNB GetVoxSizeX()
Get the voxel&#39;s size along the X axis, in mm.
#define Cerr(MESSAGE)
int ApplyOutputMaskImage()
Mask the outside of the provided input mask image.
void DeallocateImage()
Free memory for the main image matrices.
int LMS_SaveSensitivityImage(const string &a_pathToSensitivityImage, oDeformationManager *ap_DeformationManager)
FLTNB GetTimeBasisCoefficient(int a_timeBasisFunction, int a_timeFrame)
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.
int InitMultiModalImage(const vector< string > &a_pathToMultiModalImage)
int ApplyOutputFOVMasking()
Mask the outside of the transaxial FOV based on the m_fovOutPercent.
#define GEO_ROT_CCW
void DeallocateMultiModalImage()
Free memory for the multimodal image.
~oImageSpace()
oImageSpace destructor.
void DeallocateSensitivityImage()
Free memory for the sensitivity image matrices.
void InstantiateImage()
Allocate memory for the main image matrices.
int PROJ_SaveProjectionImage()
Save an image of the projected data (for analytic projection)
void InstantiateBackwardImageFromDynamicBasis(int a_nbBackwardImages)
int GetNbMultiModalImages()
Get the number of additional multimodal images.
int ApplyMaskToSensitivity()
Apply the mask to the sensitivity image (only for the first thread, the image must be reduced beforeh...
int GetNbTimeBasisFunctions()
Get the number of time basis functions.
void LMS_InstantiateForwardImage()
Allocate memory for the forward image matrices (for list-mode sensitivity generation) ...
void InstantiateSensitivityImage(const string &a_pathToSensitivityImage)
int GetSPECTSpecificParameters(uint16_t *ap_nbOfProjections, uint16_t *ap_nbHeads, FLTNB *ap_acquisitionZoom, uint16_t *ap_nbOfBins, FLTNB *ap_pixSizeXY, FLTNB *&ap_angles, FLTNB *&ap_CORtoDetectorDistance, int *ap_headRotDirection)
static sOutputManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
void LMS_CopyAtnToImage()
Copy the attenuation image contained in the &#39;m2p_attenuation&#39; matrix inside the m2p_image matrix...
int InitAttenuationImage(const string &a_pathToAtnImage)
int InitImage(const string &a_pathToInitialImage, FLTNB a_value)
void DeallocateVisitedVoxelsImage()
Free memory for the image matrix containing binary information regarding which 3D voxels have been vi...
void SaveDebugImage(const string &a_name)
void DeallocateMaskImage()
Free memory for the mask image.
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)
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...
void PROJ_DeallocateProjectionImage(int a_nbProjs)
void PrepareForwardImage()
Copy current image matrix in the forward-image buffer matrix.
FLTNB GetVoxSizeY()
Get the voxel&#39;s size along the Y axis, in mm.
int SaveOutputBasisCoefficientImage(int a_iteration, int a_subset=-1)
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.
bool GetCardStaticFlag()
Get the cardiac static flag that says if the reconstruction has only one cardiac gate or not...
void InstantiateVisitedVoxelsImage()
Memory allocation and initialization for the image matrix containing binary information regarding whi...
FLTNB GetRespBasisCoefficient(int a_respBasisFunction, int a_respGate)
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...
void LMS_DeallocateSensitivityImage()
Free memory for the sensitivity image matrices (for list-mode sensitivity generation) ...
Singleton class that Instantiate and initialize the scanner object.
int SaveOutputImage(int a_iteration, int a_subset=-1)
void LMS_DeallocateImage()
Free memory for the main image matrices (for list-mode sensitivity generation)
FLTNB GetCardBasisCoefficient(int a_cardBasisFunction, int a_cardGate)
int LoadInitialImage(const string &a_pathToImage)
int IntfWriteWholeDynBasisCoeffImgFile(const string &a_pathToImg, FLTNB ****a4p_ImgMatrix, oImageDimensionsAndQuantification *ap_ID, int vb)
Main function dedicated to Interfile 6D (dynamic dims + 3D ) image writing of basis function coeffici...
void DeallocateRefImagesForDeformation()
Free memory for the buffer sensitivity image required for image-based deformation. This function is called from the Deformation Manager.
void InstantiateForwardImage()
Allocate memory for the forward image matrices.
bool GetTimeStaticFlag()
Get the time static flag that says if the reconstruction has only one frame or not.
#define KEYWORD_OPTIONAL
void PROJ_InstantiateProjectionImage(int a_nbProjs, int a_nbPixels)
bool GetFlipOutX()
Get the boolean saying if the output image must be flipped along the X axis.
void ReduceBackwardImage(int a_imageIndex, int a_timeIndex, int a_respIndex, int a_cardIndex)
bool GetFlipOutY()
Get the boolean saying if the output image must be flipped along the Y axis.
int ApplyOutputFlip()
Just flip the output image.
This class is designed to manage the image-based deformation part of the reconstruction.
void InstantiateBackwardImageFromDynamicBins()
Allocate memory for the backward image matrices and initialize them.
void LMS_DeallocateAttenuationImage()
Free memory for the Attenuation image matrices (for analytical projection or list-mode sensitivity ge...
void DeallocateOutputImage()
Free memory for the Image matrices dedicated to output writing on disk.
Interfile fields. This structure contains all the Interfile keys currently managed by CASToR Decl...
int GetNbThreadsForImageComputation()
Get the number of threads used for image operations.
int GetNb2ndMotImgsForLMS()
call the eponym function from the oDynamicDataManager object
FLTNB * AllocateMiscellaneousImage()
Allocate a new miscellaneous image on m2p_miscellaneousImages and return the pointer to this image...
int GetNbThreadsForProjection()
Get the number of threads used for projections.
oImageSpace()
oImageSpace constructor. Initialize the member variables to their default values. ...
INTNB GetNbVoxX()
Get the number of voxels along the X axis.
INTNB GetNbVoxZ()
Get the number of voxels along the Z axis.
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 InstantiateRefImagesForDeformation()
Allocate memory for the buffer sensitivity image required for image-based deformation. This function is called from the Deformation Manager.
int SaveSensitivityImage(const string &a_pathToSensitivityImage)
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.
oImageDimensionsAndQuantification * mp_ID
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)
void LMS_InstantiateSensitivityImage()
Allocate memory for the sensitivity image matrices (for list-mode sensitivity generation) ...
void DeallocateForwardImage()
Free memory for the forward image matrices.
#define Cout(MESSAGE)
int GetNbRespBasisFunctions()
Get the number of respiratory basis functions.
int IntfKeyGetValueFromFile(const string &a_pathToHeader, const string &a_key, T *ap_return, int a_nbElts, int a_mandatoryFlag)
Look for "a_nbElts" elts in the "a_pathToHeader" interfile header matching the "a_keyword" key passed...
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.
#define KEYWORD_OPTIONAL_ERROR
void InstantiateOutputImage()
Instanciate Image matrices dedicated to output writing on disk.
void LMS_DeallocateForwardImage()
Free memory for the forward image matrices (for list-mode sensitivity generation) ...