CASToR  3.2
Tomographic Reconstruction (PET/SPECT/CT)
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 int oImageSpace::InitAttenuationImage(const string& a_pathToAtnImage)
1111 {
1112  // todo : read the number of dynamic images from Interfile header
1113  // Allocate memory and initialize matrices according to the number of fr/rg/cg
1114 
1115  // Allocate only if an image has been provided
1116  if (a_pathToAtnImage.empty()) return 0;
1117  // Verbose
1118  if (m_verbose>=VERBOSE_NORMAL) Cout("oImageSpace::InitAttenuationImage() -> Initialize and read attenuation image from file '" << a_pathToAtnImage << "'"<< endl);
1119  // Allocate only if not already done
1120  if (!m4p_attenuation)
1121  {
1123  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1124  {
1126  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
1127  {
1128  m4p_attenuation[fr][rg] = new FLTNB*[mp_ID->GetNb2ndMotImgsForLMS()];
1129 
1130  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
1131  m4p_attenuation[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
1132  }
1133  }
1134  }
1135 
1136  // Open file
1137  ifstream image_file(a_pathToAtnImage.c_str(), ios::in|ios::binary); // Read the corresponding attenuation image
1138  // Check opening
1139  if (!image_file.is_open())
1140  {
1141  Cerr("***** oImageSpace::InitAttenuationImage() -> Failed to open attenuation image from file '" << a_pathToAtnImage << "' !" << endl);
1142  return 1;
1143  }
1144  else
1145  {
1146  // Interfile image reading (INTF_LERP_ENABLED = interpolation allowed)
1147  if (IntfReadImage(a_pathToAtnImage, m4p_attenuation, mp_ID, m_verbose, INTF_LERP_ENABLED) )
1148  {
1149  Cerr("***** oImageSpace::InitAttenuationImage() -> An error occurred while reading from file '" << a_pathToAtnImage << "' !" << endl);
1150  return 1;
1151  }
1152 
1153  // Check if we have to copy the attenuation images to potential other frames/gates
1154 
1155  // Recover the number of dynamic images in the attenuation imag file
1156  int nb_atn_frames = 1, nb_atn_rgates = 1, nb_atn_cgates = 1;
1157 
1158  if( (IntfKeyGetValueFromFile(a_pathToAtnImage, "number of time frames := " , &nb_atn_frames, 1, KEYWORD_OPTIONAL) == KEYWORD_OPTIONAL_ERROR )
1159  ||(IntfKeyGetValueFromFile(a_pathToAtnImage, "number of respiratory gates := ", &nb_atn_rgates, 1, KEYWORD_OPTIONAL) == KEYWORD_OPTIONAL_ERROR )
1160  ||(IntfKeyGetValueFromFile(a_pathToAtnImage, "number of cardiac gates := " , &nb_atn_cgates, 1, KEYWORD_OPTIONAL) == KEYWORD_OPTIONAL_ERROR ) )
1161  {
1162  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);
1163  return 1;
1164  }
1165 
1166 
1167  // If more than one frame, Copy attenuation image to each other frames
1168 
1169  if(nb_atn_cgates == 1 && mp_ID->GetNb2ndMotImgsForLMS()>1)
1170  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1171  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
1172  for(int cg=1 ; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
1173  for(int v=0 ; v<mp_ID->GetNbVoxXYZ(); v++)
1174  m4p_attenuation[fr][rg][cg][v] = m4p_attenuation[fr][rg][0][v];
1175 
1176 
1177 
1178  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1179  if(nb_atn_rgates == 1 && mp_ID->GetNb1stMotImgsForLMS(fr)>1)
1180  for(int rg=1 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
1181  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
1182  for(int v=0 ; v<mp_ID->GetNbVoxXYZ(); v++)
1183  m4p_attenuation[fr][rg][cg][v] = m4p_attenuation[fr][0][cg][v];
1184 
1185 
1186 
1187  if(nb_atn_frames == 1 && mp_ID->GetNbTimeFrames()>1)
1188  for(int fr=1 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1189  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
1190  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
1191  for(int v=0 ; v<mp_ID->GetNbVoxXYZ(); v++)
1192  // If chronological motion correction is enabled, the number of motion image may be different for each frame
1193  // Check that here, just copy the attenuation image of the first motion image if it is the case
1194  m4p_attenuation[fr][rg][cg][v] = ( mp_ID->GetNb1stMotImgsForLMS(fr) == mp_ID->GetNb1stMotImgsForLMS(0) ) ?
1195  m4p_attenuation[0][rg][cg][v] :
1196  m4p_attenuation[0][0][cg][v] ;
1197 
1198 
1199  }
1200 
1201  // Close file
1202  image_file.close();
1203  return 0;
1204 }
1205 
1206 // =====================================================================
1207 // ---------------------------------------------------------------------
1208 // ---------------------------------------------------------------------
1209 // =====================================================================
1210 
1212 {
1213  // Verbose
1214  if (m_verbose>=VERBOSE_NORMAL) Cout("oImageSpace::InitAttenuationImage() -> Initialize attenuation image to uniform value " << a_value << endl);
1215  // Allocate only if not already done
1216  if (!m4p_attenuation)
1217  {
1219  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1220  {
1221  m4p_attenuation[fr] = new FLTNB**[mp_ID->GetNbRespGates()];
1222  for (int rg=0 ; rg<mp_ID->GetNbRespGates(); rg++)
1223  {
1224  m4p_attenuation[fr][rg] = new FLTNB*[mp_ID->GetNbCardGates()];
1225  for(int cg=0 ; cg<mp_ID->GetNbCardGates(); cg++)
1226  m4p_attenuation[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
1227  }
1228  }
1229  }
1230  // Initialize
1231  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1232  for (int rg=0 ; rg<mp_ID->GetNbRespGates() ; rg++)
1233  for (int cg=0 ; cg<mp_ID->GetNbCardGates() ; cg++)
1234  for (INTNB v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1235  m4p_attenuation[fr][rg][cg][v] = a_value;
1236 
1237  // Return
1238  return 0;
1239 }
1240 
1241 // =====================================================================
1242 // ---------------------------------------------------------------------
1243 // ---------------------------------------------------------------------
1244 // =====================================================================
1245 /*
1246  \fn InitImage
1247  \param a_pathToInitialImage : path to an existing image
1248  \param a_value : value to initialize each voxel with, if an input image is not provided
1249  \brief Initialize the main image, either using:
1250  - an existing image at path 'a_pathToInitialImage'
1251  - initialize each voxel with 'a_value'.
1252  \return 0 if success, positive value otherwise
1253 */
1254 int oImageSpace::InitImage(const string& a_pathToInitialImage, FLTNB a_value)
1255 {
1256  if(m_verbose>=3) Cout("oImageSpace::InitImage ..."<< endl);
1257 
1258  if (!a_pathToInitialImage.empty()) // Image initiale has been provided
1259  {
1260  if (LoadInitialImage(a_pathToInitialImage) )
1261  {
1262  Cerr("***** oImageSpace::InitImage()-> Error while trying to read file at " << a_pathToInitialImage << endl);
1263  return 1;
1264  }
1265  }
1266 
1267  else // Uniform initialization
1268  {
1269  // We initialize the content of the cylindrical FOV with the provided
1270  // value, and the exterior with a ratio of the provided value
1271  FLTNB exterior_fov_value = a_value * 1.; // For the moment we let the same value everywhere
1272  // Precast half the number of voxels over X and Y minus 1 (for efficiency)
1273  FLTNB flt_base_x = 0.5*((FLTNB)(mp_ID->GetNbVoxX()-1));
1274  FLTNB flt_base_y = 0.5*((FLTNB)(mp_ID->GetNbVoxY()-1));
1275  // Compute FOV elipse radius over X and Y, then squared
1276  FLTNB squared_radius_x = 0.5 * ((FLTNB)(mp_ID->GetNbVoxX())) * mp_ID->GetVoxSizeX();
1277  squared_radius_x *= squared_radius_x;
1278  FLTNB squared_radius_y = 0.5 * ((FLTNB)(mp_ID->GetNbVoxY())) * mp_ID->GetVoxSizeY();
1279  squared_radius_y *= squared_radius_y;
1280  // We assume that the computation of the distance from the center for a given
1281  // voxel and comparing it with the output FOV percentage costs more than performing
1282  // the loops in an inverse order compared to how the image is stored in memory.
1283  // Thus we begin the loops over X, then Y, then we test and if test passes, we
1284  // do the remaining loop over Z and over all dynamic dimensions.
1285  int x;
1286  #pragma omp parallel for private(x) schedule(guided)
1287  for (x=0; x<mp_ID->GetNbVoxX(); x++)
1288  {
1289  // Compute X distance from image center, then squared
1290  FLTNB squared_distance_x = (((FLTNB)x)-flt_base_x) * mp_ID->GetVoxSizeX();
1291  squared_distance_x *= squared_distance_x;
1292  // Loop over Y
1293  for (int y=0; y<mp_ID->GetNbVoxY(); y++)
1294  {
1295  // Compute Y distance from image center, then squared
1296  FLTNB squared_distance_y = (((FLTNB)y)-flt_base_y) * mp_ID->GetVoxSizeY();
1297  squared_distance_y *= squared_distance_y;
1298  // The value to affect to the voxel
1299  FLTNB affected_value = 0.;
1300  // Test if the voxel is inside the FOV elipse, then we affect the provided value
1301  if ( squared_distance_x/squared_radius_x + squared_distance_y/squared_radius_y <= 1. ) affected_value = a_value;
1302  // Else if outside, we affect the exterior value
1303  else affected_value = exterior_fov_value;
1304  // Loop over Z
1305  for (int z=0; z<mp_ID->GetNbVoxZ(); z++)
1306  {
1307  // Compute global voxel index
1308  INTNB index = z*mp_ID->GetNbVoxXY() + y*mp_ID->GetNbVoxX() + x;
1309  // Loops on dynamic dimensions
1310  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1311  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1312  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1313  // Affect image value
1314  m4p_image[tbf][rbf][cbf][index] = affected_value;
1315  }
1316  }
1317  }
1318  }
1319 
1320  return 0;
1321 }
1322 
1323 
1324 
1325 
1326 // =====================================================================
1327 // ---------------------------------------------------------------------
1328 // ---------------------------------------------------------------------
1329 // =====================================================================
1330 /*
1331  \fn LoadInitialImage()
1332  \param a_pathToImage : path to an existing image
1333  \brief Load the initial image provided by the user in the corresponding matrix
1334  \return 0 if success, positive value otherwise
1335 */
1336 int oImageSpace::LoadInitialImage(const string& a_pathToImage)
1337 {
1338  if(m_verbose>=3) Cout("oImageSpace::LoadInitialImage ..."<< endl);
1339 
1340  ifstream image_file;
1341  image_file.open(a_pathToImage.c_str(), ios::in | ios::binary);
1342 
1343  if(!image_file.is_open())
1344  {
1345  Cerr("***** oImageSpace::LoadInitialImage()-> Error reading file!" << endl);
1346  return 1;
1347  }
1348  else
1349  {
1350  // Interfile image reading (INTF_LERP_DISABLED = no interpolation allowed)
1351  // if(IntfReadImage(a_pathToImage, m4p_image, mp_ID, m_verbose, INTF_LERP_DISABLED))
1353  {
1354  Cerr("***** oImageSpace::LoadInitialImage()-> Error reading Interfile : " << a_pathToImage << " !" << endl);
1355  return 1;
1356  }
1357  }
1358  image_file.close();
1359 
1360  return 0;
1361 }
1362 
1363 
1364 
1365 
1366 // =====================================================================
1367 // ---------------------------------------------------------------------
1368 // ---------------------------------------------------------------------
1369 // =====================================================================
1370 /*
1371  \fn InitializationBackwardImage()
1372  \brief Initialize each voxel of the backward images to 0, also for sensitivity if not loaded (estimated on the fly).
1373 */
1375 {
1376  if(m_verbose>=3) Cout("oImageSpace::InitBackwardImage ..." << endl);
1377 
1378  // Reset backward images to 0.
1379  for (int img=0; img<m_nbBackwardImages; img++)
1380  {
1381  int th;
1382  #pragma omp parallel for private(th) schedule(static, 1)
1383  for (th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
1384  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1385  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1386  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1387  for (INTNB v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1388  {
1389  m6p_backwardImage[img][th][tbf][rbf][cbf][v] = 0.;
1390  }
1391  }
1392 
1393  // If on-the-fly sensitivity, then reset to 0.
1394  if (!m_loadedSensitivity)
1395  {
1396  int th;
1397  #pragma omp parallel for private(th) schedule(static, 1)
1398  for (th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
1399  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1400  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1401  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1402  for (INTNB v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1403  m5p_sensitivity[th][fr][rg][cg][v] = 0.;
1404  }
1405 }
1406 
1407 
1408 
1409 
1410 
1411 
1412 // =====================================================================
1413 // ---------------------------------------------------------------------
1414 // ---------------------------------------------------------------------
1415 // =====================================================================
1416 /*
1417  \fn InitSensitivityImage()
1418  \param a_pathToSensitivityImage : path to the sensitivity image (should be provided only in list-mode reconstruction)
1419  \brief Initialization for the sensitivity image matrices
1420  \details Sensitivity image initialization depends on the reconstruction mode :
1421  - list-mode: Sensitivity image has been computed before reconstruction, it is loaded from the path provided in parameter
1422  - histogram: Sensitivity image is calculated on the fly during reconstruction.
1423  First dimension (thread) is only used in histogram mode, as the on-the-fly sensitivity image computation must be thread safe
1424  \return 0 if success, positive value otherwise
1425 */
1426 int oImageSpace::InitSensitivityImage(const string& a_pathToSensitivityImage)
1427 {
1428  if(m_verbose>=3) Cout("oImageSpace::InitSensitivityImage ..."<< endl);
1429 
1430  // =====================================================================================================
1431  // Case 1: a path to a sensitivity image is provided, meaning that we are in listmode and do not compute
1432  // the sensitivity on-the-fly. In that case, we do not take the multi-threading into account.
1433  // =====================================================================================================
1434 
1435  if (!a_pathToSensitivityImage.empty())
1436  {
1437  // Open file
1438  ifstream input_file;
1439  input_file.open(a_pathToSensitivityImage.c_str(), ios::binary | ios::in);
1440  // Read file
1441  if (input_file.is_open())
1442  {
1443  // Interfile image reading (INTF_LERP_DISABLED = no interpolation allowed)
1444  if(IntfReadImage(a_pathToSensitivityImage, m5p_sensitivity[0], mp_ID, m_verbose, INTF_LERP_DISABLED) )
1445  {
1446  Cerr("***** oImageSpace::InitSensitivityImage()-> Error reading Interfile : " << a_pathToSensitivityImage << " !" << endl);
1447  return 1;
1448  }
1449  input_file.close();
1450  m_loadedSensitivity = true;
1451  }
1452  else
1453  {
1454  Cerr("***** oImageSpace::InitSensitivityImage() -> Input sensitivity file '" << a_pathToSensitivityImage << "' is missing or corrupted !" << endl);
1455  return 1;
1456  }
1457  }
1458 
1459  // =====================================================================================================
1460  // Case 2: no sensitivity provided, thus we are in histogram mode and will commpute it on-the-fly. We
1461  // thus need it to be thread-safe, so we allocated for all threads.
1462  // =====================================================================================================
1463 
1464  else
1465  {
1466 /*
1467  // Standard initialization
1468  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
1469  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1470  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1471  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1472  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1473  m5p_sensitivity[th][fr][rg][cg][v] = -1.;
1474 * */
1475  // No need to initialize to any value, as it will be done by the InitBackwardImage function which is called at the beginning
1476  // of each subset, and deals with the sensitivity initialization to 0 when m_loadedSensitivity is false
1477  m_loadedSensitivity = false;
1478  }
1479 
1480  // End
1481  return 0;
1482 }
1483 
1484 
1485 
1486 
1487 
1488 // =====================================================================
1489 // ---------------------------------------------------------------------
1490 // ---------------------------------------------------------------------
1491 // =====================================================================
1492 /*
1493  \fn oImageSpace::InitRefImagesForDeformation
1494  \brief Initialize the references image dedicated to image-based deformation. This function is called from the Deformation Manager
1495 */
1497 {
1498  if(m_verbose>=3) Cout("oDeformationManager::InitBwdImageForDeformation ..." << endl);
1499 
1500  // forward image
1501  for (int img=0; img<m_nbBackwardImages; img++)
1502  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1503  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1504  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1505  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1506  m4p_refDynForwardImage[tbf][rbf][cbf][v] = m4p_forwardImage[tbf][rbf][cbf][v];
1507 
1508  // backward image
1509  for (int img=0; img<m_nbBackwardImages; img++)
1510  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1511  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1512  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1513  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1514  m5p_refDynBackwardImage[img][tbf][rbf][cbf][v] = 0.;
1515 
1516  if(m_loadedSensitivity == false)
1517  {
1518  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1519  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1520  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1521  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
1522  m4p_refDynSensitivityImage[fr][rg][cg][v] = 0.;
1523  }
1524 }
1525 
1526 
1527 
1528 
1529 
1530 // =====================================================================
1531 // ---------------------------------------------------------------------
1532 // ---------------------------------------------------------------------
1533 // =====================================================================
1534 /*
1535  \fn InitializationSensImageForDeformation()
1536  \brief Initialize the buffer sensitivity image dedicated to image-based deformation, if required (histogram mode, as sensitivity is not loaded)
1537 
1538 void oImageSpace::InitSensImageForDeformation()
1539 {
1540  if(m_verbose>=3) Cout("oImageSpace::InitSensImageForDeformation ..."<< endl);
1541 
1542  if(m_loadedSensitivity == false)
1543  {
1544  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1545  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1546  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1547  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
1548  m4p_refDynSensitivityImage[fr][rg][cg][v] = 0.;
1549  }
1550 }
1551 */
1552 
1553 
1554 
1555 
1556 // =====================================================================
1557 // ---------------------------------------------------------------------
1558 // ---------------------------------------------------------------------
1559 // =====================================================================
1560 /*
1561  \fn InstantiateOutputImage()
1562  \brief Compute output image using the m4p_image matrix and the time/respiratory/cardiac basis functions. Store the result in the m4p_outputImage matrix
1563  \details If time/respiratory/cardiac basis functions have been initialized, this function has no effect.
1564 */
1566 {
1567  if(m_verbose>=3) Cout("oImageSpace::ComputeOutputImage ..."<< endl);
1568 
1569  // First loop on frames
1570  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1571  {
1572  // Second loop on respiratory gates
1573  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1574  {
1575  // Third loop on cardiac gates
1576  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1577  {
1578  // Reset m4p_outputImage to 0
1579  int v;
1580  #pragma omp parallel for private(v) schedule(guided)
1581  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1582  m4p_outputImage[fr][rg][cg][v] = 0.;
1583  // First loop on time basis functions
1584  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1585  {
1586  // Get frame/basis coefficient and continue if null
1587  FLTNB time_basis_coef = mp_ID->GetTimeBasisCoefficient(tbf,fr);
1588  if (time_basis_coef==0.) continue;
1589  // Second loop on respiratory basis functions
1590  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1591  {
1592  // Get resp_gate/basis coefficient and continue if null
1593  FLTNB resp_basis_coef = mp_ID->GetRespBasisCoefficient(rbf,rg);
1594  if (resp_basis_coef==0.) continue;
1595  // Third loop on cardiac basis functions
1596  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1597  {
1598  // Get card_gate_basis coefficient and continue if null
1599  FLTNB card_basis_coef = mp_ID->GetCardBasisCoefficient(cbf,cg);
1600  if (card_basis_coef==0.) continue;
1601  // Compute global coefficient
1602  FLTNB global_basis_coeff = time_basis_coef * resp_basis_coef * card_basis_coef;
1603  // Loop on voxel with OpenMP (let the chunk size by default as all values are aligned in memory)
1604  #pragma omp parallel for private(v) schedule(guided)
1605  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1606  {
1607  // Add contribution from these basis functions
1608  m4p_outputImage[fr][rg][cg][v] += m4p_image[tbf][rbf][cbf][v] * global_basis_coeff;
1609  }
1610  }
1611  }
1612  }
1613  }
1614  }
1615  }
1616 }
1617 
1618 // =====================================================================
1619 // ---------------------------------------------------------------------
1620 // ---------------------------------------------------------------------
1621 // =====================================================================
1622 
1624 {
1625  // Note: The output image is copied in the forward image at this step
1626 
1627  // If no flip, then return
1628  if ( !mp_ID->GetFlipOutX() &&
1629  !mp_ID->GetFlipOutY() &&
1630  !mp_ID->GetFlipOutZ() ) return 0;
1631  // Verbose
1632  if (m_verbose>=2)
1633  {
1634  Cout("oImageSpace::ApplyOutputFlip() -> Flip image" << endl);
1635  if (mp_ID->GetFlipOutX()) Cout(" --> Over X" << endl);
1636  if (mp_ID->GetFlipOutY()) Cout(" --> Over Y" << endl);
1637  if (mp_ID->GetFlipOutZ()) Cout(" --> Over Z" << endl);
1638  }
1639 /*
1640  // First loop on frames
1641  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1642  {
1643  // Second loop on respiratory gates
1644  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1645  {
1646  // Third loop on cardiac gates
1647  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1648  {
1649 */
1650  // First loop on time basis functions
1651  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1652  {
1653  // Second loop on respiratory basis functions
1654  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1655  {
1656  // Third loop on cardiac basis functions
1657  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1658  {
1659  // Flip over Z
1660  if (mp_ID->GetFlipOutZ())
1661  {
1662  // Half loop over Z
1663  for (INTNB z_1=0; z_1<mp_ID->GetNbVoxZ()/2; z_1++)
1664  {
1665  // Compute opposite Z
1666  INTNB z_2 = mp_ID->GetNbVoxZ() - 1 - z_1;
1667  // For efficiency
1668  INTNB base_z1 = z_1 * mp_ID->GetNbVoxXY();
1669  INTNB base_z2 = z_2 * mp_ID->GetNbVoxXY();
1670  // Loop over Y
1671  for (INTNB y=0; y<mp_ID->GetNbVoxY(); y++)
1672  {
1673  // For efficiency
1674  INTNB base_y = y * mp_ID->GetNbVoxX();
1675  // Loop over X
1676  for (INTNB x=0; x<mp_ID->GetNbVoxX(); x++)
1677  {
1678  // Compute both indices
1679  INTNB indice_1 = base_z1 + base_y + x;
1680  INTNB indice_2 = base_z2 + base_y + x;
1681  // Switch voxels
1682 /*
1683  FLTNB buffer = m4p_outputImage[fr][rg][cg][indice_1];
1684  m4p_outputImage[fr][rg][cg][indice_1] = m4p_outputImage[fr][rg][cg][indice_2];
1685  m4p_outputImage[fr][rg][cg][indice_2] = buffer;
1686 */
1687  FLTNB buffer = m4p_forwardImage[tbf][rbf][cbf][indice_1];
1688  m4p_forwardImage[tbf][rbf][cbf][indice_1] = m4p_forwardImage[tbf][rbf][cbf][indice_2];
1689  m4p_forwardImage[tbf][rbf][cbf][indice_2] = buffer;
1690  }
1691  }
1692  }
1693  }
1694  // Flip over Y
1695  if (mp_ID->GetFlipOutY())
1696  {
1697  // Loop over Z
1698  for (INTNB z=0; z<mp_ID->GetNbVoxZ(); z++)
1699  {
1700  // For efficiency
1701  INTNB base_z = z * mp_ID->GetNbVoxXY();
1702  // Half loop over Y
1703  for (INTNB y_1=0; y_1<mp_ID->GetNbVoxY()/2; y_1++)
1704  {
1705  // Compute opposite Y
1706  INTNB y_2 = mp_ID->GetNbVoxY() - 1 - y_1;
1707  // For efficiency
1708  INTNB base_y1 = y_1 * mp_ID->GetNbVoxX();
1709  INTNB base_y2 = y_2 * mp_ID->GetNbVoxX();
1710  // Loop over X
1711  for (INTNB x=0; x<mp_ID->GetNbVoxX(); x++)
1712  {
1713  // Compute both indices
1714  INTNB indice_1 = base_z + base_y1 + x;
1715  INTNB indice_2 = base_z + base_y2 + x;
1716  // Switch voxels
1717 /*
1718  FLTNB buffer = m4p_outputImage[fr][rg][cg][indice_1];
1719  m4p_outputImage[fr][rg][cg][indice_1] = m4p_outputImage[fr][rg][cg][indice_2];
1720  m4p_outputImage[fr][rg][cg][indice_2] = buffer;
1721 */
1722  FLTNB buffer = m4p_forwardImage[tbf][rbf][cbf][indice_1];
1723  m4p_forwardImage[tbf][rbf][cbf][indice_1] = m4p_forwardImage[tbf][rbf][cbf][indice_2];
1724  m4p_forwardImage[tbf][rbf][cbf][indice_2] = buffer;
1725  }
1726  }
1727  }
1728  }
1729  // Flip over X
1730  if (mp_ID->GetFlipOutX())
1731  {
1732  // Loop over Z
1733  for (INTNB z=0; z<mp_ID->GetNbVoxZ(); z++)
1734  {
1735  // For efficiency
1736  INTNB base_z = z * mp_ID->GetNbVoxXY();
1737  // Loop over Y
1738  for (INTNB y=0; y<mp_ID->GetNbVoxY(); y++)
1739  {
1740  // For efficiency
1741  INTNB base_y = y * mp_ID->GetNbVoxX();
1742  // Half loop over X
1743  for (INTNB x_1=0; x_1<mp_ID->GetNbVoxX()/2; x_1++)
1744  {
1745  // Compute opposite X
1746  INTNB x_2 = mp_ID->GetNbVoxX() - 1 - x_1;
1747  // Compute both indices
1748  INTNB indice_1 = base_z + base_y + x_1;
1749  INTNB indice_2 = base_z + base_y + x_2;
1750  // Switch voxels
1751 /*
1752  FLTNB buffer = m4p_outputImage[fr][rg][cg][indice_1];
1753  m4p_outputImage[fr][rg][cg][indice_1] = m4p_outputImage[fr][rg][cg][indice_2];
1754  m4p_outputImage[fr][rg][cg][indice_2] = buffer;
1755 */
1756  FLTNB buffer = m4p_forwardImage[tbf][rbf][cbf][indice_1];
1757  m4p_forwardImage[tbf][rbf][cbf][indice_1] = m4p_forwardImage[tbf][rbf][cbf][indice_2];
1758  m4p_forwardImage[tbf][rbf][cbf][indice_2] = buffer;
1759  }
1760  }
1761  }
1762  }
1763  }
1764  }
1765  }
1766 
1767  // Normal end
1768  return 0;
1769 }
1770 
1771 // =====================================================================
1772 // ---------------------------------------------------------------------
1773 // ---------------------------------------------------------------------
1774 // =====================================================================
1775 
1777 {
1778  // Note: The output image is copied in the forward image at this step
1779 
1780  // 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
1781  if ( mp_ID->GetFOVOutPercent()<=0. &&
1782  mp_ID->GetNbSliceOutMask()==0 ) return 0;
1783  // Verbose
1784  if (m_verbose>=2)
1785  {
1786  Cout("oImageSpace::ApplyOutputFOVMasking() -> Mask output image" << endl);
1787  if (mp_ID->GetFOVOutPercent()>0.) Cout(" --> Mask transaxial FOV outside " << mp_ID->GetFOVOutPercent() << " %" << endl);
1788  if (mp_ID->GetNbSliceOutMask()>0) Cout(" --> Mask " << mp_ID->GetNbSliceOutMask() << " slices from both axial FOV limits" << endl);
1789  }
1790  // -----------------------------------------------
1791  // Transaxial FOV masking
1792  // -----------------------------------------------
1793  if (mp_ID->GetFOVOutPercent()>0.)
1794  {
1795  // Precast half the number of voxels over X and Y minus 1 (for efficiency)
1796  FLTNB flt_base_x = 0.5*((FLTNB)(mp_ID->GetNbVoxX()-1));
1797  FLTNB flt_base_y = 0.5*((FLTNB)(mp_ID->GetNbVoxY()-1));
1798  // Compute FOV elipse radius over X and Y, then squared
1799  FLTNB squared_radius_x = 0.5 * ((FLTNB)(mp_ID->GetNbVoxX())) * mp_ID->GetVoxSizeX()
1800  * mp_ID->GetFOVOutPercent() / 100.;
1801  squared_radius_x *= squared_radius_x;
1802  FLTNB squared_radius_y = 0.5 * ((FLTNB)(mp_ID->GetNbVoxY())) * mp_ID->GetVoxSizeY()
1803  * mp_ID->GetFOVOutPercent() / 100.;
1804  squared_radius_y *= squared_radius_y;
1805  // We assume that the computation of the distance from the center for a given
1806  // voxel and comparing it with the output FOV percentage costs more than performing
1807  // the loops in an inverse order compared to how the image is stored in memory.
1808  // Thus we begin the loops over X, then Y, then we test and if test passes, we
1809  // do the remaining loop over Z and over all dynamic dimensions.
1810  int x;
1811  #pragma omp parallel for private(x) schedule(guided)
1812  for (x=0; x<mp_ID->GetNbVoxX(); x++)
1813  {
1814  // Compute X distance from image center, then squared
1815  FLTNB squared_distance_x = (((FLTNB)x)-flt_base_x) * mp_ID->GetVoxSizeX();
1816  squared_distance_x *= squared_distance_x;
1817  // Loop over Y
1818  for (int y=0; y<mp_ID->GetNbVoxY(); y++)
1819  {
1820  // Compute Y distance from image center, then squared
1821  FLTNB squared_distance_y = (((FLTNB)y)-flt_base_y) * mp_ID->GetVoxSizeY();
1822  squared_distance_y *= squared_distance_y;
1823  // Test if the voxel is inside the FOV elipse, then we skip this voxel
1824  if ( squared_distance_x/squared_radius_x + squared_distance_y/squared_radius_y <= 1. ) continue;
1825  // Loop over Z
1826  for (int z=0; z<mp_ID->GetNbVoxZ(); z++)
1827  {
1828  // Compute global voxel index
1829  INTNB index = z*mp_ID->GetNbVoxXY() + y*mp_ID->GetNbVoxX() + x;
1830 /*
1831  // First loop on frames
1832  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1833  // Second loop on respiratory gates
1834  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1835  // Third loop on cardiac gates
1836  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1837 */
1838  // First loop on time basis functions
1839  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1840  // Second loop on respiratory basis functions
1841  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1842  // Third loop on cardiac basis functions
1843  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1844  // Put 0 into this voxel
1845 // m4p_outputImage[fr][rg][cg][index] = 0.;
1846  m4p_forwardImage[tbf][rbf][cbf][index] = 0.;
1847  }
1848  }
1849  }
1850  }
1851  // -----------------------------------------------
1852  // Axial FOV masking
1853  // -----------------------------------------------
1854  if (mp_ID->GetNbSliceOutMask()>0)
1855  {
1856  INTNB removed_slices = mp_ID->GetNbSliceOutMask();
1857 /*
1858  // First loop on frames
1859  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1860  {
1861  // Second loop on respiratory gates
1862  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1863  {
1864  // Third loop on cardiac gates
1865  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1866  {
1867 */
1868  // First loop on time basis functions
1869  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1870  {
1871  // Second loop on respiratory basis functions
1872  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1873  {
1874  // Third loop on cardiac basis functions
1875  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1876  {
1877  // Mask slices
1878  for (int z=0; z<removed_slices; z++)
1879  {
1880  // First slices
1881  INTNB base_z_first = z*mp_ID->GetNbVoxXY();
1882  // Loop over Y and X
1883  for (int i=0; i<mp_ID->GetNbVoxXY(); i++)
1884  {
1885  INTNB index = base_z_first + i;
1886 // m4p_outputImage[fr][rg][cg][index] = 0.;
1887  m4p_forwardImage[tbf][rbf][cbf][index] = 0.;
1888  }
1889  // Last slices
1890  INTNB base_z_last = (mp_ID->GetNbVoxZ()-1-z)*mp_ID->GetNbVoxXY();
1891  // Loop over Y and X
1892  for (int i=0; i<mp_ID->GetNbVoxXY(); i++)
1893  {
1894  INTNB index = base_z_last + i;
1895 // m4p_outputImage[fr][rg][cg][index] = 0.;
1896  m4p_forwardImage[tbf][rbf][cbf][index] = 0.;
1897  }
1898  }
1899  }
1900  }
1901  }
1902  }
1903  // End
1904  return 0;
1905 }
1906 
1907 // =====================================================================
1908 // ---------------------------------------------------------------------
1909 // ---------------------------------------------------------------------
1910 // =====================================================================
1911 
1913 {
1914  // If no mask image return
1915  if ( !IsLoadedMask() ) return 0;
1916  // Verbose
1917  if (m_verbose>=2)
1918  {
1919  Cout("oImageSpace::ApplyOutputMaskImage() -> Mask output image with the provided input mask image" << endl);
1920  }
1921  // -----------------------------------------------
1922  // Masking
1923  // -----------------------------------------------
1924  int v;
1925  #pragma omp parallel for private(v) schedule(guided)
1926  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1927  {
1928  // First loop on frames
1929  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1930  // Second loop on respiratory gates
1931  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1932  // Third loop on cardiac gates
1933  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1934  {
1935  if (!(mp_maskImage[v]>0.)) m4p_outputImage[fr][rg][cg][v] = 0.;
1936  }
1937  }
1938 
1939  // End
1940  return 0;
1941 }
1942 
1943 // =====================================================================
1944 // ---------------------------------------------------------------------
1945 // ---------------------------------------------------------------------
1946 // =====================================================================
1947 
1949 {
1950  // If no mask image return
1951  if ( !IsLoadedMask() ) return 0;
1952 
1953  // Verbose
1954  if (m_verbose>=2)
1955  {
1956  Cout("oImageSpace::ApplyMaskToSensitivity() -> Mask the sensitivity image with the provided input mask image" << endl);
1957  }
1958  // -----------------------------------------------
1959  // Masking
1960  // -----------------------------------------------
1961  int thread_0 = 0;
1962  int v;
1963  #pragma omp parallel for private(v) schedule(guided)
1964  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1965  {
1966  // First loop on frames
1967  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1968  // Second loop on respiratory gates
1969  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1970  // Third loop on cardiac gates
1971  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1972  {
1973  if (!(mp_maskImage[v]>0.)) m5p_sensitivity[thread_0][fr][rg][cg][v] = 0.;
1974  }
1975  }
1976 
1977  // End
1978  return 0;
1979 }
1980 
1981 // =====================================================================
1982 // ---------------------------------------------------------------------
1983 // ---------------------------------------------------------------------
1984 // =====================================================================
1985 
1986 int oImageSpace::ApplyMaskToBackwardImage(int a_imageIndex, int a_timeIndex, int a_respIndex, int a_cardIndex)
1987 {
1988  // If no mask image return
1989  if ( !IsLoadedMask() ) return 0;
1990 
1991  // Verbose
1992  if (m_verbose>=2)
1993  {
1994  Cout("oImageSpace::ApplyMaskToBackwardImage() -> Mask the backward image with the provided input mask image" << endl);
1995  }
1996  // -----------------------------------------------
1997  // Masking
1998  // -----------------------------------------------
1999  int thread_0 = 0;
2000  int v;
2001  #pragma omp parallel for private(v) schedule(guided)
2002  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2003  {
2004  // First loop on frames
2005  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2006  // Second loop on respiratory gates
2007  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
2008  // Third loop on cardiac gates
2009  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
2010  {
2011  if (!(mp_maskImage[v]>0.)) m6p_backwardImage[a_imageIndex][thread_0][a_timeIndex][a_respIndex][a_cardIndex][v] = 0.;
2012  }
2013  }
2014 
2015  // End
2016  return 0;
2017 }
2018 
2019 // =====================================================================
2020 // ---------------------------------------------------------------------
2021 // ---------------------------------------------------------------------
2022 // =====================================================================
2023 /*
2024  \fn SaveOutputImage
2025  \param a_iteration : current iteration index
2026  \param a_subset : current number of subsets (or -1 by default)
2027  \brief Call the interfile function to write output image on disk
2028  \return 0 if success, positive value otherwise
2029 */
2030 int oImageSpace::SaveOutputImage(int a_iteration, int a_subset)
2031 {
2032  if (m_verbose>=3) Cout("oImageSpace::SaveOutputImage ..."<< endl);
2033 
2034  // Get the output manager
2035  sOutputManager* p_output_manager = sOutputManager::GetInstance();
2036 
2037  // Get base name
2038  string path_to_img = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
2039 
2040  // Add a suffix for iteration
2041  if (a_iteration >= 0)
2042  {
2043  stringstream ss; ss << a_iteration + 1;
2044  path_to_img.append("_it").append(ss.str());
2045  }
2046 
2047  // Add a suffix for subset (if not negative by default), this means that we save a 'subset' image
2048  if (a_subset >= 0)
2049  {
2050  stringstream ss; ss << a_subset + 1;
2051  path_to_img.append("_ss").append(ss.str());
2052  }
2053 
2054  // We need one "mode" parameter which indicates if we would like to write 1 image file or several
2055  // Adjust functions according to that
2056  if (IntfWriteImgFile(path_to_img, m4p_outputImage, mp_ID, m_verbose) )
2057  {
2058  Cerr("***** oImageSpace::SaveOutputImage()-> Error writing Interfile of output image !" << endl);
2059  return 1;
2060  }
2061 
2062  // Normal end
2063  return 0;
2064 }
2065 
2066 // =====================================================================
2067 // ---------------------------------------------------------------------
2068 // ---------------------------------------------------------------------
2069 // =====================================================================
2070 
2071 int oImageSpace::SaveOutputBasisCoefficientImage(int a_iteration, int a_subset)
2072 {
2073  if (m_verbose>=3) Cout("oImageSpace::SaveOutputBasisCoefficientImage ..."<< endl);
2074 
2075  // Get the output manager
2076  sOutputManager* p_output_manager = sOutputManager::GetInstance();
2077 
2078  // Get base name
2079  string path_to_img = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
2080 
2081  // Add a suffix for iteration
2082  if (a_iteration >= 0)
2083  {
2084  stringstream ss; ss << a_iteration + 1;
2085  path_to_img.append("_it").append(ss.str());
2086  }
2087 
2088  // Add a suffix for subset (if not negative by default), this means that we save a 'subset' image
2089  if (a_subset >= 0)
2090  {
2091  stringstream ss; ss << a_subset + 1;
2092  path_to_img.append("_ss").append(ss.str());
2093  }
2094 
2095  // At this step, the m4p_forwardImage contains the actual image basis coefficients for output
2097  {
2098  Cerr("***** oImageSpace::SaveOutputImage()-> Error writing Interfile of output image !" << endl);
2099  return 1;
2100  }
2101 
2102  // Normal end
2103  return 0;
2104 }
2105 
2106 
2107 
2108 // =====================================================================
2109 // ---------------------------------------------------------------------
2110 // ---------------------------------------------------------------------
2111 // =====================================================================
2112 /*
2113  \fn SaveDebugImage
2114  \param a_name : output name of the image
2115  \brief Just a debug function dedicated to write any kind of image on disk in raw format, for debugging purposes
2116 */
2117 void oImageSpace::SaveDebugImage(const string& a_name)
2118 {
2119  #ifdef CASTOR_DEBUG
2120  if(m_verbose>=3) Cout("oImageSpace::SaveDebugImage ..."<< endl);
2121  #endif
2122 
2123  ofstream output_file;
2124  output_file.open(a_name.c_str(), ios::binary | ios::out);
2125 
2126  //for (int fr=0; fr<ap_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
2127  // for (int rg=0; rg<ap_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
2128  // for (int cg=0; cg<ap_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
2129  // {
2130  // Write file
2131  output_file.write(reinterpret_cast<char*>(m6p_backwardImage[0][0][0][0][0]), mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
2132  // Close file
2133  output_file.close();
2134  // }
2135  // }
2136  //}
2137 }
2138 
2139 
2140 
2141 
2142 // =====================================================================
2143 // ---------------------------------------------------------------------
2144 // ---------------------------------------------------------------------
2145 // =====================================================================
2146 /*
2147  \fn PrepareForwardImage
2148  \brief Copy current image matrix in the forward-image buffer matrix
2149 */
2151 {
2152  if(m_verbose>=3) Cout("oImageSpace::PrepareForwardImage ..."<< endl);
2153 
2154  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2155  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2156  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2157  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2158  m4p_forwardImage[tbf][rbf][cbf][v] = m4p_image[tbf][rbf][cbf][v];
2159 }
2160 
2161 
2162 
2163 
2164 // =====================================================================
2165 // ---------------------------------------------------------------------
2166 // ---------------------------------------------------------------------
2167 // =====================================================================
2168 /*
2169  \fn Reduce
2170  \brief Merge parallel results into the matrix of the backward image matrix of the first thread. Also for MPI.
2171  \todo Choose computation alternative (do not know yet which one is faster...)
2172 */
2173 void oImageSpace::Reduce()
2174 {
2175  #if defined(CASTOR_OMP) || defined(CASTOR_MPI)
2176  if (m_verbose>=3) Cout("oImageSpace::Reduce() -> Merge parallel results" << endl);
2177  #endif
2178 
2179  // --------------------------------------------------------------------------------
2180  // Step 1: merge multi-threads results
2181  // --------------------------------------------------------------------------------
2182 
2183  #ifdef CASTOR_OMP
2184  if (m_verbose>=3) Cout(" --> Over threads ..." << endl);
2185 
2186  // Special case here where it appears to be always beneficial to use as many threads
2187  // as the number of threads used for projections. So we set it here, and set it back
2188  // at the end.
2189  omp_set_num_threads(mp_ID->GetNbThreadsForProjection());
2190 
2191  // todo : Choose computation alternative (do not know yet which one is faster...)
2192  int alternative = 2;
2193 
2194  // Alternative 1: standard loops based on the pointers hierarchy but mono-thread
2195  if (alternative==1 || mp_ID->GetNbThreadsForImageComputation()==1)
2196  {
2197  for (int img=0; img<m_nbBackwardImages; img++)
2198  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
2199  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2200  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2201  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2202  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2203  m6p_backwardImage[img][0][tbf][rbf][cbf][v] += m6p_backwardImage[img][th][tbf][rbf][cbf][v];
2204  }
2205  // Alternative 2: multi-thread merging with the voxels loop at first to be thread safe
2206  // Maybe faster than alternative 1, but the first loop on voxels breaks the memory order... have to try
2207  else if (alternative==2)
2208  {
2209  int v;
2210  #pragma omp parallel for private(v) schedule(guided)
2211  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2212  {
2213  for (int img=0; img<m_nbBackwardImages; img++)
2214  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
2215  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2216  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2217  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2218  m6p_backwardImage[img][0][tbf][rbf][cbf][v] += m6p_backwardImage[img][th][tbf][rbf][cbf][v];
2219  }
2220  }
2221  // If on-the-fly sensitivity, then do it also for it.
2222  if (!m_loadedSensitivity)
2223  {
2224  if (alternative==1 || mp_ID->GetNbThreadsForImageComputation()==1)
2225  {
2226  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
2227  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2228  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
2229  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
2230  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2231  m5p_sensitivity[0][fr][rg][cg][v] += m5p_sensitivity[th][fr][rg][cg][v];
2232  }
2233  else if (alternative==2)
2234  {
2235  int v;
2236  #pragma omp parallel for private(v) schedule(guided)
2237  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2238  {
2239  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
2240  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2241  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
2242  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
2243  m5p_sensitivity[0][fr][rg][cg][v] += m5p_sensitivity[th][fr][rg][cg][v];
2244  }
2245  }
2246  }
2247  // Set back the number of threads to the one for image computation
2248  omp_set_num_threads(mp_ID->GetNbThreadsForImageComputation());
2249  #endif
2250 
2251  // --------------------------------------------------------------------------------
2252  // Step 2: merge multi-instance results (MPI)
2253  // --------------------------------------------------------------------------------
2254 
2255  #ifdef CASTOR_MPI
2256  // We use the MPI_IN_PLACE as the send buffer so that it is the same as the result buffer
2257  if (m_verbose>=3) Cout(" --> Over instances ..." << endl);
2258 
2259  // Merge backward images
2260  for (int img=0; img<m_nbBackwardImages; img++)
2261  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2262  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2263  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2264  MPI_Allreduce( MPI_IN_PLACE, &m6p_backwardImage[img][0][tbf][rbf][cbf][0], mp_ID->GetNbVoxXYZ(), FLTNBMPI, MPI_SUM, MPI_COMM_WORLD );
2265 
2266  // If on-the-fly sensitivity, then do it also for it
2267  if (!m_loadedSensitivity)
2268  {
2269  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2270  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
2271  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
2272  MPI_Allreduce( MPI_IN_PLACE, &m5p_sensitivity[0][fr][rg][cg][0], mp_ID->GetNbVoxXYZ(), FLTNBMPI, MPI_SUM, MPI_COMM_WORLD );
2273  }
2274  #endif
2275 }
2276 
2277 
2278 
2279 
2280 
2281 // =====================================================================
2282 // ---------------------------------------------------------------------
2283 // ---------------------------------------------------------------------
2284 // =====================================================================
2285 /*
2286  \fn CleanNeverVisitedVoxels
2287  \brief Based on the visitedVoxelsImage, clean the never visited voxels in the image.
2288  This function must be called at the end of each iteration
2289 */
2291 {
2292  if(m_verbose>=3) Cout("oImageSpace::CleanNeverVisitedVoxels ..."<< endl);
2293  // Multi-threaded loop over voxels
2294  INTNB v;
2295  #pragma omp parallel for private(v) schedule(guided)
2296  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
2297  {
2298  // Clean this image voxel if the voxel was never visited
2299  if (mp_visitedVoxelsImage[v]==0.)
2300  {
2301  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
2302  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
2303  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
2304  m4p_image[tbf][rbf][cbf][v] = 0.;
2305  }
2306  // Reset the image
2307  mp_visitedVoxelsImage[v] = 0.;
2308  }
2309 }
2310 
2311 
2312 
2313 
2314 
2315 
2316 
2317 
2318 
2319 // LIST-MODE SENSITIVITY GENERATION FUNCTIONS
2320 // =====================================================================
2321 // ---------------------------------------------------------------------
2322 // ---------------------------------------------------------------------
2323 // =====================================================================
2324 /*
2325  \fn LMS_InstantiateImage()
2326  \brief Allocate 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_InstantiateImage ..."<< endl);
2332 
2333  m4p_image = new FLTNB***[mp_ID->GetNbTimeFrames()];
2334 
2335  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2336  {
2337  m4p_image[fr] = new FLTNB**[mp_ID->GetNb1stMotImgsForLMS(fr)];
2338 
2339  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2340  {
2341  m4p_image[fr][rg] = new FLTNB*[mp_ID->GetNb2ndMotImgsForLMS()];
2342 
2343  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2344  m4p_image[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
2345  }
2346  }
2347 }
2348 
2349 
2350 
2351 
2352 // =====================================================================
2353 // ---------------------------------------------------------------------
2354 // ---------------------------------------------------------------------
2355 // =====================================================================
2356 /*
2357  \fn LMS_DeallocateImage()
2358  \brief Free memory for the main image matrices (for list-mode sensitivity generation)
2359  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2360 */
2362 {
2363  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateImage ..."<< endl);
2364 
2365  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2366  {
2367  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2368  {
2369  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2370  if(m4p_image[fr][rg][cg] != NULL) delete m4p_image[fr][rg][cg];
2371 
2372  if(m4p_image[fr][rg] != NULL) delete[] m4p_image[fr][rg];
2373  }
2374  if(m4p_image[fr] != NULL) delete[] m4p_image[fr] ;
2375  }
2376 
2377  if(m4p_image != NULL) delete[] m4p_image;
2378  m4p_image = NULL;
2379 }
2380 
2381 
2382 
2383 // =====================================================================
2384 // ---------------------------------------------------------------------
2385 // ---------------------------------------------------------------------
2386 // =====================================================================
2387 /*
2388  \fn LMS_InstantiateForwardImage()
2389  \brief Allocate memory for the forward image matrices (for list-mode sensitivity generation)
2390  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2391 */
2393 {
2394  if(m_verbose>=3) Cout("oImageSpace::LMS_InstantiateForwardImage ..."<< endl);
2395 
2397 
2398  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2399  {
2401 
2402  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2403  {
2405 
2406  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2407  {
2408  m4p_forwardImage[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
2409  }
2410  }
2411 
2412  }
2413 }
2414 
2415 
2416 
2417 // =====================================================================
2418 // ---------------------------------------------------------------------
2419 // ---------------------------------------------------------------------
2420 // =====================================================================
2421 /*
2422  \fn LMS_DeallocateForwardImage()
2423  \brief Free memory for the forward 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_DeallocateForwardImage ..."<< endl);
2429 
2430  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2431  {
2432  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2433  {
2434  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2435  {
2436  if(m4p_forwardImage[fr][rg][cg] != NULL) delete m4p_forwardImage[fr][rg][cg];
2437  }
2438 
2439  if(m4p_forwardImage[fr][rg] != NULL) delete[] m4p_forwardImage[fr][rg];
2440  }
2441 
2442  if(m4p_forwardImage[fr] != NULL) delete[] m4p_forwardImage[fr];
2443  }
2444 
2445  if(m4p_forwardImage != NULL) delete[] m4p_forwardImage;
2446  m4p_forwardImage = NULL;
2447 }
2448 
2449 // =====================================================================
2450 // ---------------------------------------------------------------------
2451 // ---------------------------------------------------------------------
2452 // =====================================================================
2453 /*
2454  \fn LMS_InstantiateSensitivityImage()
2455  \brief Allocate memory for the sensitivity image matrices (for list-mode sensitivity generation)
2456  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2457 */
2459 {
2460  if(m_verbose>=3) Cout("oImageSpace::LMS_InstantiateSensitivityImage ..."<< endl);
2461 
2462  m5p_sensitivity = new FLTNB****[1];
2463  m5p_sensitivity[0] = new FLTNB***[mp_ID->GetNbTimeFrames()];
2464  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2465  {
2466  m5p_sensitivity[0][fr] = new FLTNB**[mp_ID->GetNb1stMotImgsForLMS(fr)];
2467  for (int rg=0; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
2468  {
2469  m5p_sensitivity[0][fr][rg] = new FLTNB*[mp_ID->GetNb2ndMotImgsForLMS()];
2470  for (int cg=0; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
2471  m5p_sensitivity[0][fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
2472  }
2473  }
2474 }
2475 
2476 
2477 
2478 
2479 // =====================================================================
2480 // ---------------------------------------------------------------------
2481 // ---------------------------------------------------------------------
2482 // =====================================================================
2483 /*
2484  \fn LMS_DeallocateSensitivityImage()
2485  \brief Free memory for the sensitivity image matrices (for list-mode sensitivity generation)
2486  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2487 */
2489 {
2490  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateSensitivityImage ..."<< endl);
2491 
2492  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2493  {
2494  for (int rg=0; rg<mp_ID->GetNb1stMotImgsForLMS(fr); rg++)
2495  {
2496  for (int cg=0; cg<mp_ID->GetNb2ndMotImgsForLMS(); cg++)
2497  {
2498  if(m5p_sensitivity[0][fr][rg][cg] != NULL) delete m5p_sensitivity[0][fr][rg][cg];
2499  }
2500  if(m5p_sensitivity[0][fr][rg] != NULL)delete m5p_sensitivity[0][fr][rg];
2501  }
2502  if(m5p_sensitivity[0][fr] != NULL)delete m5p_sensitivity[0][fr];
2503  }
2504  if(m5p_sensitivity[0] != NULL) delete[] m5p_sensitivity[0];
2505 
2506  if(m5p_sensitivity != NULL) delete[] m5p_sensitivity;
2507  m5p_sensitivity = NULL;
2508 }
2509 
2510 
2511 
2512 
2513 // =====================================================================
2514 // ---------------------------------------------------------------------
2515 // ---------------------------------------------------------------------
2516 // =====================================================================
2517 /*
2518  \fn LMS_CopyAtnToImage()
2519  \brief Copy the attenuation image contained in the 'm2p_attenuation' matrix inside the m2p_image matrix (for list-mode sensitivity generation)
2520  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2521 */
2523 {
2524  if(m_verbose>=3) Cout("oImageSpace::LMS_CopyAtnToImage ..."<< endl);
2525 
2526  if(m4p_attenuation != NULL)
2527  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2528  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2529  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2530  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2531  {
2532  m4p_image[fr][rg][cg][v] = m4p_attenuation[fr][rg][cg][v];
2533  }
2534  else
2535  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2536  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2537  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2538  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2539  {
2540  m4p_image[fr][rg][cg][v] = 0;
2541  }
2542 
2543 }
2544 
2545 
2546 // =====================================================================
2547 // ---------------------------------------------------------------------
2548 // ---------------------------------------------------------------------
2549 // =====================================================================
2550 /*
2551  \fn LMS_CopyAtnToForwardImage()
2552  \param a_use1stMotion
2553  \param a_use2ndMotion
2554  \brief Copy the attenuation image contained in the 'm2p_attenuation' matrix inside the m2p_forwardImage matrix (for list-mode sensitivity generation)
2555  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2556 */
2557 void oImageSpace::LMS_CopyAtnToForwardImage(bool a_use1stMotion, bool a_use2ndMotion)
2558 {
2559  if(m_verbose>=3) Cout("oImageSpace::LMS_CopyAtnToForwardImage ..."<< endl);
2560 
2561  if(m4p_attenuation != NULL)
2562  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2563  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2564  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2565  {
2566  // If motion is enabled, copy the attenuation image in all forward image matrices
2567  int rg_atn = a_use1stMotion ? 0 : rg ;
2568  int cg_atn = a_use1stMotion ? 0 : cg ;
2569 
2570  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2571  {
2572  m4p_forwardImage[fr][rg][cg][v] = m4p_attenuation[fr][rg_atn][cg_atn][v];
2573  }
2574  }
2575  else
2576  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2577  for(int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2578  for(int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2579  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2580  {
2581  m4p_forwardImage[fr][rg][cg][v] = 0;
2582  }
2583 }
2584 
2585 // =====================================================================
2586 // ---------------------------------------------------------------------
2587 // ---------------------------------------------------------------------
2588 // =====================================================================
2589 /*
2590  \fn LMS_CopyBackwardToSensitivityImage()
2591  \brief Copy the backward image containing the result of the sensitivity back-projection, into the sensitivity matrix.
2592  \details This function is dedicated to list-mode sensitivity (LMS) generation, it is useful to apply convolution and to save the image.
2593 */
2595 {
2596  if(m_verbose>=3) Cout("oImageSpace::LMS_CopyBackwardToSensitivity ..."<< endl);
2597  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2598  for (int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2599  for (int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2600  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2601  m5p_sensitivity[0][fr][rg][cg][v] = m6p_backwardImage[0][0][fr][rg][cg][v];
2602 }
2603 
2604 
2605 // =====================================================================
2606 // ---------------------------------------------------------------------
2607 // ---------------------------------------------------------------------
2608 // =====================================================================
2609 /*
2610  \fn LMS_PrepareForwardImage()
2611  \brief Copy current image in forward-image buffer (for list-mode sensitivity generation)
2612  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2613 */
2615 {
2616  if(m_verbose>=3) Cout("oImageSpace::LMS_PrepareForwardImage ..."<< endl);
2617 
2618  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2619  for (int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2620  for (int cg=0 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2621  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2622  m4p_forwardImage[fr][rg][cg][v] = m4p_image[fr][rg][cg][v];
2623 }
2624 
2625 
2626 
2627 
2628 // =====================================================================
2629 // ---------------------------------------------------------------------
2630 // ---------------------------------------------------------------------
2631 // =====================================================================
2632 
2633 void oImageSpace::ReduceBackwardImage(int a_imageIndex, int a_timeIndex, int a_respIndex, int a_cardIndex)
2634 {
2635  #if defined(CASTOR_OMP) || defined(CASTOR_MPI)
2636  if (m_verbose>=3) Cout("oImageSpace::ReduceBackwardImage() -> Merge parallel results" << endl);
2637  #endif
2638 
2639  // --------------------------------------------------------------------------------
2640  // Step 1: merge multi-threads results
2641  // --------------------------------------------------------------------------------
2642 
2643  #ifdef CASTOR_OMP
2644  // Verbose
2645  if (m_verbose>=3) Cout(" --> Over threads ..." << endl);
2646  // Do it
2647  for (int th=1 ; th<mp_ID->GetNbThreadsForProjection() ; th++)
2648  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2649  m6p_backwardImage[a_imageIndex][0][a_timeIndex][a_respIndex][a_cardIndex][v] += m6p_backwardImage[a_imageIndex][th][a_timeIndex][a_respIndex][a_cardIndex][v];
2650  #endif
2651 
2652  // --------------------------------------------------------------------------------
2653  // Step 2: merge multi-instance results (MPI)
2654  // --------------------------------------------------------------------------------
2655 
2656  #ifdef CASTOR_MPI
2657  // We use the MPI_IN_PLACE as the send buffer so that it is the same as the result buffer
2658  if (m_verbose>=3) Cout(" --> Over instances ..." << endl);
2659  // Merge backward images
2660  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 );
2661  #endif
2662 }
2663 
2664 
2665 
2666 
2667 // =====================================================================
2668 // ---------------------------------------------------------------------
2669 // ---------------------------------------------------------------------
2670 // =====================================================================
2671 /*
2672  \fn LMS_SaveSensitivityImage
2673  \param a_pathToSensitivityImage : path to the sensitivity image (should be provided only in list-mode reconstruction)
2674  \param ap_DeformationManager : Pointer to the deformation manager objet (required to retrieve the number of gates in the sensitivity image)
2675  \brief Call the interfile function to write the sensitivity image on disk
2676  \details If image deformation is enabled for respiratory/cardiac gated data, the gated images are summed up into one image and normalize
2677  \return 0 if success, positive value otherwise
2678  \todo Interfile management : if the number of sensitivity image to write is different to the actual number of image
2679  as it could be the case for dynamic imaging, if one sensitivity image is required for the whole dynamic serie,
2680  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)
2681 */
2682 int oImageSpace::LMS_SaveSensitivityImage(const string& a_pathToSensitivityImage, oDeformationManager* ap_DeformationManager)
2683 {
2684  if(m_verbose>=3) Cout("oImageSpace::LMS_SaveSensitivityImage ..."<< endl);
2685 
2686  // If image-based motion is enabled, list-mode sensitivity images code generates a set of several provisional images of the reference position.
2687  // Each image is generated using one different forward/backward deformations parameters which will be used during reconstruction.
2688  // - For gated (respiratory/cardiac) motion correction, a simple average of these images is performed here for each frame.
2689  // - For timestamp-based motion (involuntary patient motion), the averaging depends on which transformation occurs in each frame.
2690  // The provisional sensitivity images are weighs accordingly to the duration of each motion subset inside the frame.
2691  // (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)
2692 
2693  int nb_reco_card_images = ap_DeformationManager->GetNbSensImagesCardDeformation(mp_ID->GetNb2ndMotImgsForLMS());
2694 
2695  // Average of the cardiac images, if cardiac gating is enabled (no entry in the loop of cardiac images otherwise)
2696  if (ap_DeformationManager->GetNbSensImagesCardDeformation(mp_ID->GetNb2ndMotImgsForLMS()) == 0)
2697  {
2698  nb_reco_card_images = 1;
2699  FLTNB card_normalization = mp_ID->GetNb2ndMotImgsForLMS();
2700 
2701  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2702  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames(); fr++)
2703  for (int rg=0 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2704  {
2705  for (int cg=1 ; cg<mp_ID->GetNb2ndMotImgsForLMS() ; cg++)
2706  m5p_sensitivity[0][fr][rg][0][v] += m5p_sensitivity[0][fr][rg][cg][v];
2707  m5p_sensitivity[0][fr][rg][0][v] /= card_normalization;
2708  }
2709  }
2710 
2711  // Average of the respiratory images, if respiratory gating is enabled (no entry in the loop of respiratory images)
2712  if (!mp_ID->IsPMotionEnabled() // No patient motion correction
2713  && ap_DeformationManager->GetNbSensImagesRespDeformation(mp_ID->GetNb1stMotImgsForLMS(0)) == 0)
2714  {
2715  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames(); fr++)
2716  {
2717  FLTNB resp_normalization = mp_ID->GetNb1stMotImgsForLMS(fr);
2718 
2719  for (int cg=0 ; cg<nb_reco_card_images ; cg++) // Dual gating : Be sure to get the right number of card images
2720  {
2721  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2722  {
2723  for (int rg=1 ; rg<mp_ID->GetNb1stMotImgsForLMS(fr) ; rg++)
2724  m5p_sensitivity[0][fr][0][cg][v] += m5p_sensitivity[0][fr][rg][cg][v];
2725 
2726  m5p_sensitivity[0][fr][0][cg][v] /= resp_normalization;
2727  }
2728  }
2729  }
2730  }
2731 
2732  // Average of the sensitivity images when involuntary patient motion correction is enabled
2733  if (mp_ID->IsPMotionEnabled() )
2734  {
2735  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames(); fr++)
2736  {
2737  //FLTNB ipm_normalization = GetNb1stMotImgsForLMS(fr);
2738 
2739  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2740  {
2741 
2742 
2743  //for (int ipm=1 ; ipm<mp_ID->GetNb1stMotImgsForLMS(fr) ; ipm++)
2744  // m5p_sensitivity[0][fr][0][0][v] += m5p_sensitivity[0][fr][ipm][0][v] ;
2745 
2746  //m5p_sensitivity[0][fr][0][0][v] /= ipm_normalization;
2747 
2748  FLTNB sensitivity_value_avg = 0.;
2749 
2750  // Add the contribution of the sensitivity images with appropriate weighing
2751  for (int ipm=0 ; ipm<mp_ID->GetNb1stMotImgsForLMS(fr) ; ipm++)
2752  sensitivity_value_avg += m5p_sensitivity[0][fr][ipm][0][v] * mp_ID->GetListPMotionWeightInFrameForLMS(fr, ipm);
2753 
2754  // Recover the average value in the first sensitivity image
2755  m5p_sensitivity[0][fr][0][0][v] = sensitivity_value_avg;
2756  }
2757 
2758  }
2759  }
2760 
2761 
2762 
2763  if (SaveSensitivityImage(a_pathToSensitivityImage) )
2764  {
2765  Cerr("***** oImageSpace::LMS_SaveSensitivityImage()-> Error writing Sensitivity image !" << endl);
2766  return 1;
2767  }
2768 
2769  return 0;
2770 }
2771 
2772 
2773 
2774 
2775 // =====================================================================
2776 // ---------------------------------------------------------------------
2777 // ---------------------------------------------------------------------
2778 // =====================================================================
2779 /*
2780  \fn SaveSensitivityImage
2781  \param a_pathToSensitivityImage : path to the sensitivity image (should be provided only in list-mode reconstruction)
2782  \brief Call the interfile function to write the sensitivity image on disk
2783  \return 0 if success, positive value otherwise
2784 */
2785 int oImageSpace::SaveSensitivityImage(const string& a_pathToSensitivityImage)
2786 {
2787  if(m_verbose>=3) Cout("oImageSpace::SaveSensitivityImage ..."<< endl);
2788 
2789  if (IntfWriteImgFile(a_pathToSensitivityImage, m5p_sensitivity[0], mp_ID, m_verbose) )
2790  {
2791  Cerr("***** oImageSpace::SaveSensitivityImage()-> Error writing Interfile of output image !" << endl);
2792  return 1;
2793  }
2794 
2795  return 0;
2796 }
2797 
2798 
2799 
2800 
2801 // PROJECTION SCRIPT FUNCTIONS
2802 // =====================================================================
2803 // ---------------------------------------------------------------------
2804 // ---------------------------------------------------------------------
2805 // =====================================================================
2806 /*
2807  \fn PROJ_InstantiateProjectionImage()
2808  \param a_nbProjs : a number of projection slices in the projection
2809  \param a_nbPixels : a total number of pixels in the projection slices
2810  \brief Instanciate and initialize projection image for analytical projection
2811  \details This function is currently only dedicated to SPECT projection, and used by the analytical projection script
2812  \todo support for PET (requires projection format)
2813 */
2814 void oImageSpace::PROJ_InstantiateProjectionImage(int a_nbProjs, int a_nbPixels)
2815 {
2816  if(m_verbose>=3) Cout("oImageSpace::PROJ_InstantiateProjectionImage ..."<< endl);
2817 
2818  m2p_projectionImage = new FLTNB*[a_nbProjs];
2819 
2820  for(int p=0 ; p<a_nbProjs ; p++)
2821  {
2822  m2p_projectionImage[p] = new FLTNB[a_nbPixels];
2823  for(int px=0 ; px<a_nbPixels ; px++)
2824  m2p_projectionImage[p][px] = 0;
2825  }
2826 }
2827 
2828 
2829 
2830 
2831 // =====================================================================
2832 // ---------------------------------------------------------------------
2833 // ---------------------------------------------------------------------
2834 // =====================================================================
2835 /*
2836  \fn PROJ_DeallocateProjectionImage()
2837  \param a_nbProjs : a number of projection slices in the projection
2838  \brief Free memory for the projection image for analytical projection
2839  \details This function is currently only dedicated to SPECT projection, and used by the analytical projection script
2840  \todo support for PET (requires projection format)
2841 */
2843 {
2844  if(m_verbose>=3) Cout("oImageSpace::PROJ_DeallocateProjectionImage ..."<< endl);
2845 
2846  for(int p=0 ; p<a_nbProjs ; p++)
2847  {
2848  delete[] m2p_projectionImage[p];
2849  }
2850 
2851  delete[] m2p_projectionImage;
2852  m2p_projectionImage = NULL;
2853 }
2854 
2855 
2856 
2857 
2858 // =====================================================================
2859 // ---------------------------------------------------------------------
2860 // ---------------------------------------------------------------------
2861 // =====================================================================
2862 /*
2863  \fn PROJ_InitImage()
2864  \param a_pathToInitialImage : path to the image to project
2865  \brief Load the initial image for the analytical projection
2866  \return 0 if success, positive value otherwise
2867 */
2868 int oImageSpace::PROJ_InitImage(const string& a_pathToInitialImage)
2869 {
2870  if(m_verbose>=3) Cout("oImageSpace::PROJ_InitImage ..."<< endl);
2871 
2872  if (!a_pathToInitialImage.empty()) // Image initiale
2873  {
2874  if(PROJ_LoadInitialImage(a_pathToInitialImage) )
2875  {
2876  Cerr("***** oImageSpace::PROJ_InitImage()-> Error while trying to read file at " << a_pathToInitialImage << endl);
2877  return 1;
2878  }
2879  }
2880  else
2881  {
2882  {
2883  Cerr("***** oImageSpace::PROJ_InitImage()-> No projected image provided ! " << endl);
2884  return 1;
2885  }
2886  }
2887  return 0;
2888 }
2889 
2890 
2891 
2892 
2893 // =====================================================================
2894 // ---------------------------------------------------------------------
2895 // ---------------------------------------------------------------------
2896 // =====================================================================
2897 /*
2898  \fn PROJ_LoadInitialImage()
2899  \param a_pathToImage : path to the image to project
2900  \brief Load the initial image for the analytical projection
2901  \return 0 if success, positive value otherwise
2902 */
2903 int oImageSpace::PROJ_LoadInitialImage(const string& a_pathToImage)
2904 {
2905  if(m_verbose>=3) Cout("oImageSpace::PROJ_LoadInitialImage ..."<< endl);
2906 
2907  ifstream image_file;
2908  image_file.open(a_pathToImage.c_str(), ios::in | ios::binary);
2909 
2910  if(!image_file.is_open())
2911  {
2912  Cerr("***** oImageSpace::PROJ_LoadInitialImage()-> Error reading file !" << endl);
2913  return 1;
2914  }
2915  else
2916  {
2917  // Interfile image reading (INTF_LERP_DISABLED = no interpolation allowed)
2919  {
2920  Cerr("***** oImageSpace::PROJ_LoadInitialImage()-> Error reading Interfile : " << a_pathToImage << " !" << endl);
2921  return 1;
2922  }
2923 
2924  }
2925  image_file.close();
2926  return 0;
2927 }
2928 
2929 
2930 
2931 
2932 
2933 // =====================================================================
2934 // ---------------------------------------------------------------------
2935 // ---------------------------------------------------------------------
2936 // =====================================================================
2937 /*
2938  \fn PROJ_SaveProjectionImage
2939  \brief Save an image of the projected data
2940  \return 0 if success, positive value otherwise
2941  \todo Support for PET projeted data
2942  implement interfile file recovery in a scanner class function
2943 */
2945 {
2946  if(m_verbose>=3) Cout("oImageSpace::PROJ_SaveProjectionImage ..."<< endl);
2947 
2948  string img_name = "_ProjectionImage";
2949 
2950  // Initialize interfile fields structure
2951  Intf_fields Img_fields;
2952  IntfKeySetFieldsOutput(&Img_fields , mp_ID);
2953 
2954  Img_fields.process_status = "acquired";
2955 
2956  // Get specific info about projection
2958 
2959  // Recover the values in the interfile structure.
2960  // todo : This step for PET projection
2961  // todo : Implement this step in scanner class functions
2962  // todo : For spect, deal with systems with several detector heads (currently assume one head)
2963  if(p_ScanMgr->GetScannerType() == 2)
2964  {
2965  uint16_t nb_projs, nb_heads;
2966  FLTNB zoom;
2967  uint16_t nb_bins[2];
2968  FLTNB pix_size[2];
2969  FLTNB* p_angles;
2970  FLTNB* p_radius;
2971  int dir_rot;
2972  p_ScanMgr->GetSPECTSpecificParameters(&nb_projs,
2973  &nb_heads,
2974  &zoom,
2975  nb_bins,
2976  pix_size,
2977  p_angles,
2978  p_radius,
2979  &dir_rot);
2980 
2981  Img_fields.nb_projections = nb_projs;
2982  Img_fields.nb_detector_heads = nb_heads;
2983  Img_fields.mtx_size[0] = nb_bins[0];
2984  Img_fields.mtx_size[1] = nb_bins[1];
2985  Img_fields.vox_size[0] = pix_size[0];
2986  Img_fields.vox_size[1] = pix_size[1];
2987  Img_fields.extent_rotation = 360; // 360 by default
2988  // Check rotation direction from the two first angles
2989  Img_fields.direction_rotation = (dir_rot == GEO_ROT_CCW)?
2990  "CCW" :
2991  "CW";
2992  Img_fields.first_angle = p_angles[0];
2993 
2994  // Note: In Interfile v3.3, doesn't seem to have a field to provide
2995  // individually each angle and Center of rotation.
2996  // Looks like it was planned for ulterior version, at least for the
2997  // center of rotation
2998  // For now, we just write each projection angle and radius as an array key
2999 
3000  string angles_str = "{";
3001  string radius_str = "{";
3002 
3003  // Flag to check if the radius is similar for each projection angle
3004  bool has_single_radius = true;
3005 
3006  for(uint16_t p=0 ; p<nb_projs ; p++)
3007  {
3008  stringstream ss_a, ss_r;
3009  // projection angles
3010  ss_a << p_angles[p];
3011  angles_str.append(ss_a.str());
3012  (p<nb_projs-1) ? angles_str.append(",") : angles_str.append("}");
3013 
3014  // radius
3015  ss_r << p_radius[p];
3016  radius_str.append(ss_r.str());
3017  (p<nb_projs-1) ? radius_str.append(",") : radius_str.append("}");
3018  if(p_radius[p] != p_radius[0])
3019  has_single_radius = false;
3020 
3021  // Some interfile editors struggle with long line, perhaps because
3022  // they are stored in char[256] or something like that
3023  // Hence, break line after a certain number of elements
3024  if((p+1)%30 == 0)
3025  {
3026  angles_str.append("\n");
3027  radius_str.append("\n");
3028  }
3029  }
3030 
3031  // If there is one common radius, write its value in the related string
3032  if(has_single_radius)
3033  {
3034  stringstream ss;
3035  ss << p_radius[0];
3036  radius_str = ss.str();
3037  }
3038 
3039  Img_fields.projection_angles = angles_str;
3040  Img_fields.radius = radius_str;
3041 
3042 
3043  // Common code to all modalities
3044  sOutputManager* p_output_manager = sOutputManager::GetInstance();
3045  string img_file_name = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
3046  img_file_name.append("_ProjImage");
3047 
3048  if(IntfWriteProjFile(img_file_name, m2p_projectionImage, mp_ID, Img_fields, m_verbose) )
3049  {
3050  Cerr("***** oImageSpace::PROJ_SaveProjectionImage()-> Error writing Interfile of output image !" << endl);
3051  return 1;
3052  }
3053  }
3054  return 0;
3055 }
3056 
3057 
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) ...