CASToR  2.0
Tomographic Reconstruction (PET/SPECT/CT)
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
iDataFilePET.cc
Go to the documentation of this file.
1 /*
2 This file is part of CASToR.
3 
4  CASToR is free software: you can redistribute it and/or modify it under the
5  terms of the GNU General Public License as published by the Free Software
6  Foundation, either version 3 of the License, or (at your option) any later
7  version.
8 
9  CASToR is distributed in the hope that it will be useful, but WITHOUT ANY
10  WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
12  details.
13 
14  You should have received a copy of the GNU General Public License along with
15  CASToR (in file GNU_GPL.TXT). If not, see <http://www.gnu.org/licenses/>.
16 
17 Copyright 2017-2018 all CASToR contributors listed below:
18 
19  --> current contributors: Thibaut MERLIN, Simon STUTE, Didier BENOIT, Claude COMTAT, Marina FILIPOVIC, Mael MILLARDET
20  --> past contributors: Valentin VIELZEUF
21 
22 This is CASToR version 2.0.
23 */
24 
32 #include "iDataFilePET.hh"
33 
34 // =====================================================================
35 // ---------------------------------------------------------------------
36 // ---------------------------------------------------------------------
37 // =====================================================================
38 
40 {
41  // Set all members to default values
45  m_maxAxialDiffmm = -1.;
46  m_isotope = "unknown";
47  m_resolutionTOF = -1.;
48  m_nbBinsTOF = -1;
49  m_binSizeTOF = -1.;
50  m_eventKindFlag = false;
51  m_atnCorrectionFlag = false;
53  m_normCorrectionFlag = false;
55  m_scatCorrectionFlag = false;
57  m_randCorrectionFlag = false;
59  m_TOFInfoFlag = false;
60  m_ignoreTOFFlag = false;
62 }
63 
64 // =====================================================================
65 // ---------------------------------------------------------------------
66 // ---------------------------------------------------------------------
67 // =====================================================================
68 
70 
71 // =====================================================================
72 // ---------------------------------------------------------------------
73 // ---------------------------------------------------------------------
74 // =====================================================================
75 
76 int iDataFilePET::ReadSpecificInfoInHeader(bool a_affectQuantificationFlag)
77 {
78  // Verbose
80  if (m_verbose>=VERBOSE_DETAIL) Cout("iDataFilePET::ReadSpecificInfoInHeader() -> Read information specific to PET" << endl);
81  // Read optional fields in the header, check if errors (issue during data reading/conversion (==1) )
82  if (ReadDataASCIIFile(m_headerFileName, "Maximum number of lines per event", &m_maxNumberOfLinesPerEvent, 1, KEYWORD_OPTIONAL)==1 ||
83  ReadDataASCIIFile(m_headerFileName, "Maximum axial difference mm", &m_maxAxialDiffmm, 1, KEYWORD_OPTIONAL)==1 ||
84  ReadDataASCIIFile(m_headerFileName, "TOF resolution (ps)", &m_resolutionTOF, 1, KEYWORD_OPTIONAL)==1 ||
85  ReadDataASCIIFile(m_headerFileName, "TOF information flag", &m_TOFInfoFlag, 1, KEYWORD_OPTIONAL)==1 ||
86  ReadDataASCIIFile(m_headerFileName, "List TOF measurement range (ps)", &m_TOFMeasurementRange, 1, KEYWORD_OPTIONAL)==1 ||
87  ReadDataASCIIFile(m_headerFileName, "Histo TOF number of bins", &m_nbBinsTOF, 1, KEYWORD_OPTIONAL)==1 ||
88  ReadDataASCIIFile(m_headerFileName, "Histo TOF bin size (ps)", &m_binSizeTOF, 1, KEYWORD_OPTIONAL)==1 ||
89  ReadDataASCIIFile(m_headerFileName, "Coincidence kind flag", &m_eventKindFlag, 1, KEYWORD_OPTIONAL)==1 ||
90  ReadDataASCIIFile(m_headerFileName, "Attenuation correction flag", &m_atnCorrectionFlag, 1, KEYWORD_OPTIONAL)==1 ||
91  ReadDataASCIIFile(m_headerFileName, "Normalization correction flag", &m_normCorrectionFlag, 1, KEYWORD_OPTIONAL)==1 ||
92  ReadDataASCIIFile(m_headerFileName, "Scatter correction flag", &m_scatCorrectionFlag, 1, KEYWORD_OPTIONAL)==1 ||
93  ReadDataASCIIFile(m_headerFileName, "Random correction flag", &m_randCorrectionFlag, 1, KEYWORD_OPTIONAL)==1 ||
95  {
96  Cerr("***** iDataFilePET::ReadSpecificInfoInHeader() -> Error while reading optional fields in the header data file !" << endl);
97  return 1;
98  }
99  // Give the PET isotope to the oImageDimensionsAndQuantification that manages the quantification factors
100  if (a_affectQuantificationFlag && mp_ID->SetPETIsotope(m_bedIndex, m_isotope))
101  {
102  Cerr("***** iDataFilePET::ReadSpecificInfoInHeader() -> A problem occured while setting the isotope to oImageDimensionsAndQuantification !" << endl);
103  return 1;
104  }
105  // End
106  return 0;
107 }
108 
109 // =====================================================================
110 // ---------------------------------------------------------------------
111 // ---------------------------------------------------------------------
112 // =====================================================================
113 
115 {
116  iDataFilePET* p_DataFilePET = (dynamic_cast<iDataFilePET*>(ap_DataFile));
118  m_maxAxialDiffmm = p_DataFilePET->GetMaxAxialDiffmm();
119  m_isotope = p_DataFilePET->GetIsotope();
120  m_resolutionTOF = p_DataFilePET->GetTOFResolution();
121  m_nbBinsTOF = p_DataFilePET->GetNbTOFBins();
122  m_binSizeTOF = p_DataFilePET->GetTOFBinSize();
123  m_TOFInfoFlag = p_DataFilePET->GetTOFInfoFlag();
124  m_eventKindFlag = p_DataFilePET->GetEventKindFlag();
125  m_atnCorrectionFlag = p_DataFilePET->GetAtnCorrectionFlag();
126  m_normCorrectionFlag = p_DataFilePET->GetNormCorrectionFlag();
127  m_scatCorrectionFlag = p_DataFilePET->GetScatCorrectionFlag();
128  m_randCorrectionFlag = p_DataFilePET->GetRandCorrectionFlag();
130  // End
131  return 0;
132 }
133 
134 // =====================================================================
135 // ---------------------------------------------------------------------
136 // ---------------------------------------------------------------------
137 // =====================================================================
138 
140 {
141  // Verbose
143  if (m_verbose>=VERBOSE_DETAIL) Cout("iDataFilePET::ComputeSizeEvent() -> In bytes" << endl);
144 
145  // For MODE_LIST events
146  if (m_dataMode == MODE_LIST)
147  {
148  // Size of the mandatory element in a list-mode event: Time + max_nb_lines*2*crystalID
149  m_sizeEvent = sizeof(uint32_t) + m_maxNumberOfLinesPerEvent*2*sizeof(uint32_t);
150  // Number of lines
151  if (m_maxNumberOfLinesPerEvent>1) m_sizeEvent += sizeof(uint16_t);
152  // Optional flags
153  if (m_eventKindFlag) m_sizeEvent += sizeof(uint8_t);
154  if (m_atnCorrectionFlag) m_sizeEvent += sizeof(FLTNBDATA);
158  if (m_TOFInfoFlag) m_sizeEvent += sizeof(FLTNBDATA);
159  // POI availability in each direction
160  for (int i=0 ; i<3; i++) if (mp_POIDirectionFlag[i]) m_sizeEvent += 2*sizeof(FLTNBDATA);
161  }
162  // For MODE_HISTOGRAM events
163  else if (m_dataMode == MODE_HISTOGRAM)
164  {
165  // Size of the mandatory element in a histo event: Time + nbBinsTOF*event_value + max_nb_line*2*crystalID
166  m_sizeEvent = sizeof(uint32_t) + m_nbBinsTOF*sizeof(FLTNBDATA) + m_maxNumberOfLinesPerEvent*2*sizeof(uint32_t);
167  // Number of lines
168  if (m_maxNumberOfLinesPerEvent>1) m_sizeEvent += sizeof(uint16_t);
169  // Optional flags
170  if (m_atnCorrectionFlag) m_sizeEvent += sizeof(FLTNBDATA);
174  }
175  // For MODE_NORMALIZATION events
176  else if (m_dataMode == MODE_NORMALIZATION)
177  {
178  // max_nb_lines*2*crystalID
179  m_sizeEvent = m_maxNumberOfLinesPerEvent*2*sizeof(uint32_t);
180  // Number of lines
181  if (m_maxNumberOfLinesPerEvent>1) m_sizeEvent += sizeof(uint16_t);
182  // Optional flag for normalization
184  // Optional flag for attenuation
185  if (m_atnCorrectionFlag) m_sizeEvent += sizeof(FLTNBDATA);
186  }
187  // Unknown event type
188  else
189  {
190  Cerr("***** iDataFilePET::ComputeSizeEvent() -> Unknown event mode !" << endl);
191  return 1;
192  }
193 
194  // Check
195  if (m_sizeEvent<=0)
196  {
197  Cerr("***** iDataFilePET::ComputeSizeEvent() -> Error, the Event size in bytes should be >= 0 !" << endl;);
198  return 1;
199  }
200 
201  // Verbose
202  if (m_verbose>=VERBOSE_DETAIL) Cout(" --> Event size = " << m_sizeEvent << " bytes" << endl);
203  // End
204  return 0;
205 }
206 
207 // =====================================================================
208 // ---------------------------------------------------------------------
209 // ---------------------------------------------------------------------
210 // =====================================================================
211 
213 {
214  // Verbose
217  {
218  if (m_dataMode==MODE_HISTOGRAM) Cout("iDataFilePET::PrepareDataFile() -> Build histogram events" << endl);
219  else if (m_dataMode==MODE_LIST) Cout("iDataFilePET::PrepareDataFile() -> Build listmode events" << endl);
220  else if (m_dataMode==MODE_NORMALIZATION) Cout("iDataFilePET::PrepareDataFile() -> Build normalization events" << endl);
221  }
222 
223  // ==============================================================================
224  // Allocate event buffers (one for each thread)
225  // ==============================================================================
226  if (m_verbose>=VERBOSE_DETAIL) Cout(" --> Allocating an event buffer for each thread" << endl);
227  // Instanciation of the event buffer according to the data type
229 
230  // Allocate the events per each thread
231  for (int th=0 ; th<mp_ID->GetNbThreadsForProjection() ; th++)
232  {
233 
234  // For MODE_LIST events
235  if (m_dataMode == MODE_LIST)
236  {
237  m2p_BufferEvent[th] = new iEventListPET();
238  // If TOF, convert and transmit information about the total range of TOF measurements ( conversion from delta time into delta length )
240  {
241  ((iEventListPET*)m2p_BufferEvent[th])->SetHasTOFinfo(true);
242  ((iEventListPET*)m2p_BufferEvent[th])->SetSpatialTOFMeasurementRange(m_TOFMeasurementRange*SPEED_OF_LIGHT*0.5);
243  }
244  }
245  // For MODE_HISTOGRAM events
246  else if (m_dataMode == MODE_HISTOGRAM)
247  {
248  m2p_BufferEvent[th] = new iEventHistoPET();
249  // If we ignore the TOF information, then the event buffer will have only one TOF bin;
250  // the TOF contributions will be sum up when reading the event.
251  if (m_ignoreTOFFlag) ((iEventHistoPET*)m2p_BufferEvent[th])->SetEventNbTOFBins(1);
252  else ((iEventHistoPET*)m2p_BufferEvent[th])->SetEventNbTOFBins(m_nbBinsTOF);
253  }
254  // For MODE_NORMALIZATION events
255  else if (m_dataMode == MODE_NORMALIZATION)
256  {
257  m2p_BufferEvent[th] = new iEventNorm();
258  }
259 
260  // Set the maximum number of lines per event
262 
263  // Allocate crystal IDs
264  if (m2p_BufferEvent[th]->AllocateID())
265  {
266  Cerr("*****iDataFilePET::PrepareDataFile() -> Error while trying to allocate memory for the Event object for thread " << th << " !" << endl;);
267  return 1;
268  }
269  }
270 
271  // ==============================================================================
272  // Deal with specific corrections
273  // ==============================================================================
274 
275  // In case of attenuation correction flag, see if we ignore this correction
277  {
278  // Affect the ignored flag from the ignored corrections list processed by the oImageDimensionsAndQuantification
280  // Verbose
282  {
283  if (m_ignoreAttnCorrectionFlag) Cout(" --> Ignore attenuation correction" << endl);
284  else Cout(" --> Correct for attenuation" << endl);
285  }
286  }
287  // In case of normalization correction flag, see if we ignore this correction
289  {
290  // Affect the ignored flag from the ignored corrections list processed by the oImageDimensionsAndQuantification
292  // Verbose
294  {
295  if (m_ignoreNormCorrectionFlag) Cout(" --> Ignore normalization correction" << endl);
296  else Cout(" --> Correct for normalization" << endl);
297  }
298  }
299  // In case of scatter correction flag, see if we ignore this correction
301  {
302  // Affect the ignored flag from the ignored corrections list processed by the oImageDimensionsAndQuantification
304  // Verbose
306  {
307  if (m_ignoreScatCorrectionFlag) Cout(" --> Ignore scatter correction" << endl);
308  else Cout(" --> Correct for scatter events" << endl);
309  }
310  }
311  // In case of random correction flag, see if we ignore this correction
313  {
314  // Affect the ignored flag from the ignored corrections list processed by the oImageDimensionsAndQuantification
316  // Verbose
318  {
319  if (m_ignoreRandCorrectionFlag) Cout(" --> Ignore random correction" << endl);
320  else Cout(" --> Correct for random events" << endl);
321  }
322  }
323 
324  // ==============================================================================
325  // Deal with TOF
326  // ==============================================================================
327 
328  // TOF information availability and use
329  if (m_TOFInfoFlag)
330  {
331  // Special case when TOF is ignored but we are in list-mode, we have scatter correction factors and we do not ignore them.
332  // In this case, the scatter correction will be completely wrong, so we make the code crash here
334  {
335  Cerr("***** iDataFilePET::PrepareDataFile() -> Scatter correction for list-mode TOF data will be erroneous if TOF information is ignored in the projections !" << endl);
336  return 1;
337  }
338  // Verbose
340  {
341  if (m_ignoreTOFFlag) Cout(" --> TOF information available but ignored " << endl);
342  else Cout(" --> Use TOF information" << endl);
343  }
344  }
345 
346  // Normal end
347  return 0;
348 }
349 
350 // =====================================================================
351 // ---------------------------------------------------------------------
352 // ---------------------------------------------------------------------
353 // =====================================================================
354 
355 vEvent* iDataFilePET::GetEventSpecific(char* ap_buffer, int a_th)
356 {
358 
359  // Work on a copy of the input pointer
360  char* file_position = ap_buffer;
361 
362  // For MODE_LIST PET data
363  if (m_dataMode == MODE_LIST)
364  {
365  // Cast the event pointer
366  iEventListPET* event = (dynamic_cast<iEventListPET*>(m2p_BufferEvent[a_th]));
367  // Mandatory time field: [uint32_t (time)]
368  event->SetTimeInMs(*reinterpret_cast<uint32_t*>(file_position));
369  file_position += sizeof(uint32_t);
370  // Optional attenuation correction field: [FLTNBDATA (attnCorrectionFactor)]
372  {
373  if (!m_ignoreAttnCorrectionFlag) event->SetAttenuationCorrectionFactor(*reinterpret_cast<FLTNBDATA*>(file_position));
374  file_position += sizeof(FLTNBDATA);
375  }
376  // Optional kind: [uint8_t kind]
377  if (m_eventKindFlag)
378  {
379  event->SetKind(*reinterpret_cast<uint8_t*>(file_position));
380  file_position += sizeof(uint8_t);
381  }
382  // Optional scatter correction field: [FLTNBDATA (scatter)]
384  {
385  if (!m_ignoreScatCorrectionFlag) event->SetScatterRate(0, *reinterpret_cast<FLTNBDATA*>(file_position));
386  file_position += sizeof(FLTNBDATA);
387  }
388  // Optional random correction field: [FLTNBDATA (random)]
390  {
391  if (!m_ignoreRandCorrectionFlag) event->SetRandomRate(*reinterpret_cast<FLTNBDATA*>(file_position));
392  file_position += sizeof(FLTNBDATA);
393  }
394  // Optional normalization correction field: [FLTNBDATA (norm)]
396  {
397  if (!m_ignoreNormCorrectionFlag) event->SetNormalizationFactor(*reinterpret_cast<FLTNBDATA*>(file_position));
398  file_position += sizeof(FLTNBDATA);
399  }
400  // Optional TOF information field: [FLTNBDATA (TOF)]
401  if (m_TOFInfoFlag)
402  {
403  if (!m_ignoreTOFFlag) event->SetTOFMeasurement(*reinterpret_cast<FLTNBDATA*>(file_position));
404  file_position += sizeof(FLTNBDATA);
405  }
406  // Optional POI correction fields: [FLTNBDATA (POI1[3])] [FLTNBDATA (POI2[3])]
407  for (int i=0; i<3; i++)
408  {
409  if (mp_POIDirectionFlag[i])
410  {
411  if (!m_ignorePOIFlag) event->SetPOI1(i,*reinterpret_cast<FLTNBDATA*>(file_position));
412  file_position += sizeof(FLTNBDATA);
413  if (!m_ignorePOIFlag) event->SetPOI2(i,*reinterpret_cast<FLTNBDATA*>(file_position));
414  file_position += sizeof(FLTNBDATA);
415  }
416  }
417  // Mandatory nb contributing LORs if m_maxNumberOfLinesPerEvent>1: [uint16_t nbLines]
419  {
420  event->SetNbLines(*reinterpret_cast<uint16_t*>(file_position));
421  file_position += sizeof(uint16_t);
422  }
423  // Mandatory crystal IDs: [uint32_t (ID1)] [uint32_t (ID2)]
424  for (int i=0 ; i<event->GetNbLines() ; i++)
425  {
426  event->SetID1(i, *reinterpret_cast<uint32_t*>(file_position));
427  file_position += sizeof(uint32_t);
428  event->SetID2(i, *reinterpret_cast<uint32_t*>(file_position));
429  file_position += sizeof(uint32_t);
430  }
431  }
432 
433  // For MODE_HISTOGRAM PET data
434  else if (m_dataMode == MODE_HISTOGRAM)
435  {
436  // Cast the event pointer
437  iEventHistoPET* event = (dynamic_cast<iEventHistoPET*>(m2p_BufferEvent[a_th]));
438  // Mandatory time field: [uint32_t (time)]
439  event->SetTimeInMs(*reinterpret_cast<uint32_t*>(file_position));
440  file_position += sizeof(uint32_t);
441  // Optional attenuation correction field: [FLTNBDATA (attnCorrectionFactor)]
443  {
444  if (!m_ignoreAttnCorrectionFlag) event->SetAttenuationCorrectionFactor(*reinterpret_cast<FLTNBDATA*>(file_position));
445  file_position += sizeof(FLTNBDATA);
446  }
447 
448  // Optional random correction field: [FLTNBDATA (random)]
450  {
451  if (!m_ignoreRandCorrectionFlag) event->SetRandomRate(*reinterpret_cast<FLTNBDATA*>(file_position));
452  file_position += sizeof(FLTNBDATA);
453  }
454 
455  // Optional normalization correction field: [FLTNBDATA (norm)]
457  {
458  if (!m_ignoreNormCorrectionFlag) event->SetNormalizationFactor(*reinterpret_cast<FLTNBDATA*>(file_position));
459  file_position += sizeof(FLTNBDATA);
460  }
461  // If we ignore the TOF information, then we have to sum up the bins' contributions
462  if (m_ignoreTOFFlag)
463  {
464  // Compute the total of event values and scatter rates
465  FLTNBDATA total_event_value = 0.;
466  FLTNBDATA total_scatter_rate = 0.;
467  for (int tb=0 ; tb<m_nbBinsTOF ; tb++)
468  {
469  // Mandatory bin value: [FLTNBDATA bin value]
470  total_event_value += *reinterpret_cast<FLTNBDATA*>(file_position);
471  file_position += sizeof(FLTNBDATA);
472  // Optional scatter correction field: [FLTNBDATA (scatter)]
474  {
475  total_scatter_rate += *reinterpret_cast<FLTNBDATA*>(file_position);
476  file_position += sizeof(FLTNBDATA);
477  }
478  }
479  // Affect the total to the event
480  event->SetEventValue(0,total_event_value);
481  if (!m_ignoreScatCorrectionFlag) event->SetScatterRate(0,total_scatter_rate);
482  }
483  // Otherwise we set all different TOF bin contributions
484  else
485  {
486  for (int tb=0 ; tb<m_nbBinsTOF ; tb++)
487  {
488  // Mandatory bin value: [FLTNBDATA bin value]
489  event->SetEventValue(tb, *reinterpret_cast<FLTNBDATA*>(file_position));
490  file_position += sizeof(FLTNBDATA);
491  // Optional scatter correction field: [FLTNBDATA (scatter)]
493  {
494  if (!m_ignoreScatCorrectionFlag) event->SetScatterRate(tb, *reinterpret_cast<FLTNBDATA*>(file_position));
495  file_position += sizeof(FLTNBDATA);
496  }
497  }
498  }
499  // Mandatory nb contributing LORs if m_maxNumberOfLinesPerEvent>1: [uint16_t nbLines]
501  {
502  event->SetNbLines(*reinterpret_cast<uint16_t*>(file_position));
503  file_position += sizeof(uint16_t);
504  }
505  // Mandatory crystal IDs: [uint32_t (c1)] [uint32_t (c2)]
506  for (int i=0 ; i<event->GetNbLines() ; i++)
507  {
508  event->SetID1(i, *reinterpret_cast<uint32_t*>(file_position));
509  file_position += sizeof(uint32_t);
510  event->SetID2(i, *reinterpret_cast<uint32_t*>(file_position));
511  file_position += sizeof(uint32_t);
512  }
513  }
514 
515  // For MODE_NORMALIZATION PET data
516  else if (m_dataMode == MODE_NORMALIZATION)
517  {
518  // Cast the event pointer
519  iEventNorm* event = (dynamic_cast<iEventNorm*>(m2p_BufferEvent[a_th]));
520  // Optional attenuation correction field: [FLTNBDATA (norm)]
522  {
523  if (!m_ignoreAttnCorrectionFlag) event->SetAttenuationCorrectionFactor(*reinterpret_cast<FLTNBDATA*>(file_position));
524  file_position += sizeof(FLTNBDATA);
525  }
526  // Optional normalization correction field: [FLTNBDATA (norm)]
528  {
529  if (!m_ignoreNormCorrectionFlag) event->SetNormalizationFactor(*reinterpret_cast<FLTNBDATA*>(file_position));
530  file_position += sizeof(FLTNBDATA);
531  }
532  // Mandatory nb contributing LORs if m_maxNumberOfLinesPerEvent>1: [uint16_t nbLines]
534  {
535  event->SetNbLines(*reinterpret_cast<uint16_t*>(file_position));
536  file_position += sizeof(uint16_t);
537  }
538  // Mandatory crystal IDs: [uint32_t (c1)] [uint32_t (c2)]
539  for (int i=0 ; i<event->GetNbLines() ; i++)
540  {
541  event->SetID1(i, *reinterpret_cast<uint32_t*>(file_position));
542  file_position += sizeof(uint32_t);
543  event->SetID2(i, *reinterpret_cast<uint32_t*>(file_position));
544  file_position += sizeof(uint32_t);
545  }
546  }
547 
548  // Return the updated event
549  return m2p_BufferEvent[a_th];
550 }
551 
552 // =====================================================================
553 // ---------------------------------------------------------------------
554 // ---------------------------------------------------------------------
555 // =====================================================================
556 
558 {
560  if (m_verbose==0) return;
561  // Describe the datafile
562  Cout("iDataFilePET::DescribeSpecific() -> Here is some specific content of the PET datafile" << endl);
563  Cout(" --> Isotope: " << m_isotope << endl);
564  if (m_TOFInfoFlag)
565  {
566  Cout(" --> TOF resolution: " << m_resolutionTOF << " ps" << endl);
567  if (m_dataMode==MODE_LIST)
568  {
569  Cout(" --> TOF measurement range (for list-mode data): " << m_TOFMeasurementRange << " ps" << endl);
570  }
571  else if (m_dataMode==MODE_HISTOGRAM)
572  {
573  Cout(" --> Number of TOF bins: " << m_nbBinsTOF << endl);
574  Cout(" --> TOF bin size: " << m_binSizeTOF << " ps" << endl);
575  }
576  }
577  if (m_dataMode==MODE_LIST && m_eventKindFlag) Cout(" --> Event kind is present" << endl);
578  if (m_scatCorrectionFlag) Cout(" --> Scatter correction is present" << endl);
579  if (m_randCorrectionFlag) Cout(" --> Random correction is present" << endl);
580  if (m_atnCorrectionFlag) Cout(" --> Attenuation correction is present" << endl);
581  if (m_normCorrectionFlag) Cout(" --> Normalization correction is present" << endl);
582  if (m_maxNumberOfLinesPerEvent>1) Cout(" --> Maximum number of lines per event: " << m_maxNumberOfLinesPerEvent << endl);
583 }
584 
585 // =====================================================================
586 // ---------------------------------------------------------------------
587 // ---------------------------------------------------------------------
588 // =====================================================================
589 
591 {
593  // Error if m_dataType != PET
594  if (m_dataType != TYPE_PET)
595  {
596  Cerr("***** iDataFilePET::CheckSpecificParameters() -> Data type should be PET !'" << endl);
597  return 1;
598  }
599  // Check if Maximum ring difference has been initialized, use default value (negative) otherwise.
600  if (m_maxAxialDiffmm < 0.)
601  {
602  m_maxAxialDiffmm = -1.;
603  }
604  // Some checks for POI use
605  if (m_POIInfoFlag)
606  {
607  // Not compatible with histogram data
609  {
610  Cerr("***** iDataFilePET::CheckSpecificParameters() -> POI correction flag is enabled while the data are histogrammed, no sense !" << endl);
611  return 1;
612  }
613  // For each direction (radial, tangential and axial), look at the resolution. If not negative, the info should be the datafile
614  for (int i=0; i<3; i++) if (mp_POIResolution[i]>=0.) mp_POIDirectionFlag[i] = true;
615  }
616  // Some checks for TOF use
617  if (m_TOFInfoFlag)
618  {
619  // Check if the resolution was provided
620  if (m_resolutionTOF<0.)
621  {
622  Cerr("***** iDataFilePET::CheckSpecificParameters() -> TOF information is used while there is no TOF resolution specified in the datafile header !" << endl);
623  return 1;
624  }
625  // For histogram data
627  {
628  // Check if the number of bins has been provided
629  if (m_nbBinsTOF<=1)
630  {
631  Cerr("***** iDataFilePET::CheckSpecificParameters() -> TOF information is used while there is only one TOF bin (specified in the datafile header) !" << endl);
632  return 1;
633  }
634  // Check if the bin size has been provided
635  if (m_binSizeTOF<=0.)
636  {
637  Cerr("***** iDataFilePET::CheckSpecificParameters() -> TOF information is used while there is no bin size specified in the datafile header !" << endl);
638  return 1;
639  }
640  }
641  // For list mode data
642  else if (m_dataMode==MODE_LIST)
643  {
644  // set the number of TOF bins to 1
645  m_nbBinsTOF = 1;
646  return 0;
647  }
648  }
649  else
650  {
651  // Set the number of TOF bins to 1
652  m_nbBinsTOF = 1;
653  }
654  // End
655  return 0;
656 }
657 
658 // =====================================================================
659 // ---------------------------------------------------------------------
660 // ---------------------------------------------------------------------
661 // =====================================================================
662 
664 {
666 
667  // Create the stream
668  fstream* p_file = new fstream( m_dataFileName.c_str(), ios::binary| ios::in );
669  // Check that datafile exists
670  if (!p_file->is_open())
671  {
672  Cerr("***** iDataFilePET::CheckFileSizeConsistency() -> Failed to open input file '" << m_dataFileName.c_str() << "' !" << endl);
673  Cerr(" (Provided in the data file header: " << m_headerFileName << ")" << endl);
674  return 1;
675  }
676  // Get file size in bytes
677  p_file->seekg(0, ios::end);
678  int64_t sizeInBytes = p_file->tellg();
679  // Close stream and delete it
680  p_file->close();
681  delete p_file;
682  // Check datafile self-consistency
683  if (m_nbEvents*m_sizeEvent != sizeInBytes)
684  {
685  Cerr("--------------------------------------------------------------------------------------------------------------------------------------" << endl);
686  Cerr("***** iDataFilePET::CheckFileSizeConsistency() -> DataFile size is not consistent with the information provided by the user/datafile !" << endl);
687  Cerr(" --> Expected size : "<< m_nbEvents*m_sizeEvent << endl);
688  Cerr(" --> Actual size : "<< sizeInBytes << endl << endl);
689  Cerr(" ADDITIONAL INFORMATION ABOUT THE DATAFILE INITIALIZATION : " << endl);
690  if (m_maxNumberOfLinesPerEvent > 1) Cerr(" --> Compression enabled. Each event should contain " << m_maxNumberOfLinesPerEvent <<
691  "LORs, or an equivalent number of LOR + garbage LORs" << endl);
692  else Cerr(" --> No compression in the data" << endl);
693  if (m_eventKindFlag) Cerr(" --> Coincidence kind term is enabled" << endl);
694  else Cerr(" --> No information about the kind of coincidences in the data" << endl);
695  if (m_normCorrectionFlag) Cerr(" --> Normalization correction term is enabled" << endl);
696  else Cerr(" --> No normalization correction term in the data" << endl);
697  if (m_atnCorrectionFlag) Cerr(" --> Attenuation correction term is enabled" << endl);
698  else Cerr(" --> No attenuation correction term in the data" << endl);
699  if (m_scatCorrectionFlag) Cerr(" --> Scatter correction term is enabled" << endl);
700  else Cerr(" --> No scatter correction term in the data" << endl);
701  if (m_randCorrectionFlag) Cerr(" --> Random correction term is enabled" << endl);
702  else Cerr(" --> No random correction term in the data" << endl);
703  if (m_TOFInfoFlag) Cerr(" --> TOF information is enabled. Resolution is: " << m_resolutionTOF << " ps"<< endl);
704  else Cerr(" --> No TOF information in the data" << endl);
705  Cerr(" --> Calibration factor value is: " << m_calibrationFactor << endl);
706  if (mp_POIResolution[0]<0.) Cerr(" --> No POI enabled on the radial axis" << endl);
707  else Cerr(" --> POI resolution on the radial axis is: " << mp_POIResolution[0] << endl);
708  if (mp_POIResolution[1]<0.) Cerr(" --> No POI enabled on the tangential axis" << endl);
709  else Cerr(" --> POI resolution on the tangential axis is: " << mp_POIResolution[1] << endl);
710  if (mp_POIResolution[2]<0.) Cerr(" --> No POI enabled on the axial axis" << endl);
711  else Cerr(" --> POI resolution on the axial axis is: " << mp_POIResolution[2] << endl);
712  Cerr("--------------------------------------------------------------------------------------------------------------------------------------" << endl);
713  return 1;
714  }
715  // End
716  return 0;
717 }
718 
719 // =====================================================================
720 // ---------------------------------------------------------------------
721 // ---------------------------------------------------------------------
722 // =====================================================================
723 
725 {
727  // Dynamic cast the vDataFile to a iDataFilePET
728  iDataFilePET* p_data_file = (dynamic_cast<iDataFilePET*>(ap_DataFile));
729  // Check maximum number of lines per event
731  {
732  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Maximum numbers of lines by event are inconsistent !" << endl);
733  return 1;
734  }
735  // Check maximum ring difference
736  if (m_maxAxialDiffmm!=p_data_file->GetMaxAxialDiffmm())
737  {
738  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Maximum ring differences are inconsistent !" << endl);
739  return 1;
740  }
741  // Check isotope
742  if (m_isotope!=p_data_file->GetIsotope())
743  {
744  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Isotopes are inconsistent !" << endl);
745  return 1;
746  }
747  // Check event kind flag
748  if (m_eventKindFlag!=p_data_file->GetEventKindFlag())
749  {
750  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Event kind flags are inconsistent !" << endl);
751  return 1;
752  }
753  // Check attenuation correction flag
754  if (m_atnCorrectionFlag!=p_data_file->GetAtnCorrectionFlag())
755  {
756  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Attenuation correction flags are inconsistent !" << endl);
757  return 1;
758  }
759  // Check normalization correction flag
760  if (m_normCorrectionFlag!=p_data_file->GetNormCorrectionFlag())
761  {
762  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Normalization correction flags are inconsistent !" << endl);
763  return 1;
764  }
765  // Check scatter correction flag
766  if (m_scatCorrectionFlag!=p_data_file->GetScatCorrectionFlag())
767  {
768  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Scatter correction flags are inconsistent !" << endl);
769  return 1;
770  }
771  // Check random correction flag
772  if (m_randCorrectionFlag!=p_data_file->GetRandCorrectionFlag())
773  {
774  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Random correction flags are inconsistent !" << endl);
775  return 1;
776  }
777  // Check TOF info flag
778  if (m_TOFInfoFlag!=p_data_file->GetTOFInfoFlag())
779  {
780  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> TOF information flags are inconsistent !" << endl);
781  return 1;
782  }
783  // Check TOF resolution
784  if (m_resolutionTOF!=p_data_file->GetTOFResolution())
785  {
786  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> TOF resolutions are inconsistent !" << endl);
787  return 1;
788  }
789  // Check number of TOF bins
790  if (m_nbBinsTOF!=p_data_file->GetNbTOFBins())
791  {
792  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Numbers of TOF bins are inconsistent !" << endl);
793  return 1;
794  }
795  // Check TOF bin size
796  if (m_binSizeTOF!=p_data_file->GetTOFBinSize())
797  {
798  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> TOF bin sizes are inconsistent !" << endl);
799  return 1;
800  }
801  // Check data mode
802  if (m_dataMode!=p_data_file->GetDataMode())
803  {
804  Cerr("***** iDataFilePET::CheckSpecificConsistencyWithAnotherDataFile() -> Data modes are inconsistent (list-mode or histogram) !" << endl);
805  return 1;
806  }
807  // End
808  return 0;
809 }
810 
811 // =====================================================================
812 // ---------------------------------------------------------------------
813 // ---------------------------------------------------------------------
814 // =====================================================================
815 
817 {
819  if (m_verbose>=VERBOSE_DETAIL) Cout("iDataFilePET::PROJ_InitFile() -> Initialize datafile for writing"<< endl);
820 
821  m_nbEvents = 0; //Std initialization for projection
822  m_nbBinsTOF = 1; // Initialization of this variable required for PROJ_Write functions
823 
824  m_isotope = "unknown";
825 
826  // Default time initialization if image dimensions object has not been initialized
827  if(mp_ID->IsInitialized())
828  {
830 
831  for(uint16_t fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
833 
834  // Instanciate a fstream datafile for each thread
835  m2p_dataFile = new fstream*[mp_ID->GetNbThreadsForProjection()];
836  }
837  else
838  {
839  m_startTimeInSec = 0.;
840  m_durationInSec = 1.;
841  m2p_dataFile = new fstream*[1];
842  }
843 
844 
845 
846  // Set name of the projection data header
847  sOutputManager* p_OutMgr;
848  p_OutMgr = sOutputManager::GetInstance();
849  string path_name = p_OutMgr->GetPathName();
850  string img_name = p_OutMgr->GetBaseName();
851  m_headerFileName = path_name.append(img_name).append("_CstrProj").append(".Cdh");
852 
853  for (int th=0 ; th<mp_ID->GetNbThreadsForProjection() ; th++)
854  {
855  m_dataFileName = m_headerFileName.substr(0, m_headerFileName.find_last_of(".")).append(".Cdf"); // Generate datafile name from header file
856 
857  // We'll write projeted data in several files (corresponding to the number of thread) to be concatenated at the end of the projection process.
858  stringstream ss;
859  ss << th;
860  string datafile_name = m_dataFileName;
861  datafile_name.append("_").append(ss.str());
862 
863  m2p_dataFile[th] = new fstream( datafile_name.c_str(), ios::binary | ios::out | ios::trunc);
864  }
865 
866 
867  // Check there is no existing datafile with such name. Remove it otherwise
868  ifstream fcheck(m_dataFileName.c_str());
869  if(fcheck.good())
870  {
871  fcheck.close();
872  #ifdef _WIN32
873  string dos_instruction = "del " + m_dataFileName;
874  system(dos_instruction.c_str());
875  #else
876  remove(m_dataFileName.c_str());
877  #endif
878  }
879 
880  if (m_verbose>=VERBOSE_DETAIL) Cout(" --> Output datafile name :" << m_dataFileName << endl);
881 
882  return 0;
883 }
884 
885 // =====================================================================
886 // ---------------------------------------------------------------------
887 // ---------------------------------------------------------------------
888 // =====================================================================
889 
890 int iDataFilePET::WriteEvent(vEvent* ap_Event, int a_th)
891 {
893 
894  if (m_dataMode == MODE_LIST)
895  {
896  if (WriteListEvent((iEventListPET*)ap_Event, a_th))
897  {
898  Cerr("*****iDataFilePET::WriteEvent() -> Error while trying to write projection datafile (list-mode)" << endl;);
899  return 1;
900  }
901  }
902 
903  if (m_dataMode == MODE_HISTOGRAM)
904  {
905  if (WriteHistoEvent((iEventHistoPET*)ap_Event, a_th))
906  {
907  Cerr("*****iDataFilePET::WriteEvent() -> Error while trying to write projection datafile (histogram)" << endl;);
908  return 1;
909  }
910  }
911 
912  return 0;
913 }
914 
915 // =====================================================================
916 // ---------------------------------------------------------------------
917 // ---------------------------------------------------------------------
918 // =====================================================================
919 
921 {
923 
924  // Write sequentially each field of the event according to the type of the event.
925  m2p_dataFile[a_th]->clear();
926 
927  uint32_t time = ap_Event->GetTimeInMs();
928  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&time), sizeof(uint32_t));
929 
931  {
932  FLTNB atn_corr_factor = ap_Event->GetAtnCorrFactor();
933  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&atn_corr_factor), sizeof(FLTNB));
934  }
935 
937  {
938  FLTNB rdm_rate = ap_Event->GetEventRdmRate();
939  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&rdm_rate), sizeof(FLTNB));
940  }
941 
943  {
944  FLTNB norm_corr_factor = ap_Event->GetNormFactor();
945  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&norm_corr_factor), sizeof(FLTNB));
946  }
947 
948  for(int b=0 ; b<m_nbBinsTOF ; b++)
949  {
950  FLTNB event_value = ap_Event->GetEventValue(b);
951  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&event_value), sizeof(FLTNB));
953  {
954  FLTNB scat_rate = ap_Event->GetEventScatRate(b);
955  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&scat_rate), sizeof(FLTNB));
956  }
957  }
958 
959  uint16_t nb_lines = ap_Event->GetNbLines();
960  // Write the number of lines only if the m_maxNumberOfLinesPerEvent is above 1
961  if (m_maxNumberOfLinesPerEvent>1) m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&nb_lines), sizeof(uint16_t));
962 
963  for(int i=0 ; i<nb_lines ; i++)
964  {
965  uint32_t id1 = ap_Event->GetID1(i);
966  uint32_t id2 = ap_Event->GetID2(i);
967  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&id1), sizeof(uint32_t));
968  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&id2), sizeof(uint32_t));
969  }
970 
971  // 0-filling if needed (nb lines inferior to the max number for PET data with compression)
972  if(nb_lines<m_maxNumberOfLinesPerEvent)
973  for(int i=0 ; i<m_maxNumberOfLinesPerEvent-nb_lines ; i++)
974  {
975  uint32_t gbg = 0;
976  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&gbg), sizeof(uint32_t));
977  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&gbg), sizeof(uint32_t));
978  }
979 
980  return 0;
981 }
982 
983 // =====================================================================
984 // ---------------------------------------------------------------------
985 // ---------------------------------------------------------------------
986 // =====================================================================
987 
989 {
991 
992  // Write sequentially each field of the event according to the type of the event.
993  m2p_dataFile[a_th]->clear();
994 
995  uint32_t time = ap_Event->GetTimeInMs();
996  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&time), sizeof(uint32_t));
997 
998 
1000  {
1001  FLTNB atn_corr_factor = ap_Event->GetAtnCorrFactor();
1002  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&atn_corr_factor), sizeof(FLTNB));
1003  }
1004 
1005  if(m_eventKindFlag)
1006  {
1007  uint8_t event_kind = ap_Event->GetKind();
1008  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&event_kind), sizeof(uint8_t));
1009  }
1010 
1012  {
1013  FLTNB scat_rate = ap_Event->GetEventScatRate(0);
1014  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&scat_rate), sizeof(FLTNB));
1015  }
1016 
1018  {
1019  FLTNB rdm_rate = ap_Event->GetEventRdmRate();
1020  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&rdm_rate), sizeof(FLTNB));
1021  }
1022 
1024  {
1025  FLTNB norm_corr_factor = ap_Event->GetNormFactor();
1026  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&norm_corr_factor), sizeof(FLTNB));
1027  }
1028 
1029  if(m_TOFInfoFlag)
1030  {
1031  FLTNBDATA TOF_measurement = ap_Event->GetTOFMeasurement();
1032  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&TOF_measurement), sizeof(FLTNBDATA));
1033  }
1034 
1035 
1036  if(mp_POIResolution[0]>0)
1037  {
1038  FLTNB POI_1 = ap_Event->GetPOI1(0);
1039  FLTNB POI_2 = ap_Event->GetPOI2(0);
1040  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&POI_1), sizeof(FLTNB));
1041  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&POI_2), sizeof(FLTNB));
1042  }
1043  if(mp_POIResolution[1]>0)
1044  {
1045  FLTNB POI_1 = ap_Event->GetPOI1(1);
1046  FLTNB POI_2 = ap_Event->GetPOI2(1);
1047  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&POI_1), sizeof(FLTNB));
1048  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&POI_2), sizeof(FLTNB));
1049  }
1050  if(mp_POIResolution[2]>0)
1051  {
1052  FLTNB POI_1 = ap_Event->GetPOI1(2);
1053  FLTNB POI_2 = ap_Event->GetPOI2(2);
1054  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&POI_1), sizeof(FLTNB));
1055  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&POI_2), sizeof(FLTNB));
1056  }
1057 
1058  uint16_t nb_lines = ap_Event->GetNbLines();
1059  // Write the number of lines only if the m_maxNumberOfLinesPerEvent is above 1
1060  if (m_maxNumberOfLinesPerEvent>1) m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&nb_lines), sizeof(uint16_t));
1061 
1062  for(int i=0 ; i<nb_lines ; i++)
1063  {
1064  uint32_t id1 = ap_Event->GetID1(i);
1065  uint32_t id2 = ap_Event->GetID2(i);
1066  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&id1), sizeof(uint32_t));
1067  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&id2), sizeof(uint32_t));
1068  }
1069 
1070  // 0-filling if needed (nb lines inferior to the max number for PET data with compression)
1071  if(nb_lines<m_maxNumberOfLinesPerEvent)
1072  for(int i=0 ; i<m_maxNumberOfLinesPerEvent-nb_lines ; i++)
1073  {
1074  uint32_t gbg = 0;
1075  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&gbg), sizeof(uint32_t));
1076  m2p_dataFile[a_th]->write(reinterpret_cast<char*>(&gbg), sizeof(uint32_t));
1077  }
1078 
1079  return 0;
1080 }
1081 
1082 // =====================================================================
1083 // ---------------------------------------------------------------------
1084 // ---------------------------------------------------------------------
1085 // =====================================================================
1086 
1088 {
1089  // Verbose
1091  if (m_verbose>=VERBOSE_NORMAL) Cout("iDataFilePET::WriteHeader() -> Write header file '" << m_headerFileName << "'" << endl);
1092  // Open file
1093  fstream headerFile;
1094  headerFile.open(m_headerFileName.c_str(), ios::out);
1095 
1096  if (!headerFile.is_open())
1097  {
1098  Cerr("***** iDataFilePET::WriteHeader() -> Failed to open output header file '" << m_headerFileName << "' !" << endl);
1099  return 1;
1100  }
1101 
1102  // Data file name
1103  headerFile << "Data filename: " << GetFileFromPath(m_dataFileName) << endl;
1104  // Number of events
1105  headerFile << "Number of events: " << m_nbEvents << endl;
1106  // Data mode
1107  if (m_dataMode==MODE_HISTOGRAM) headerFile << "Data mode: histogram" << endl;
1108  else if (m_dataMode==MODE_LIST) headerFile << "Data mode: list-mode" << endl;
1109  else if (m_dataMode==MODE_NORMALIZATION) headerFile << "Data mode: normalization" << endl;
1110  // PET data type
1111  headerFile << "Data type: PET" << endl;
1112  // Acquisition start time in seconds
1113  headerFile << "Start time (s): " << m_startTimeInSec << endl;
1114  // Acquisition duration in seconds
1115  headerFile << "Duration (s): " << m_durationInSec << endl;
1116  // Scanner name
1117  headerFile << "Scanner name: " << sScannerManager::GetInstance()->GetScannerName() << endl;
1118  // Maximum axial difference
1119  if(m_maxAxialDiffmm > 0.)
1120  headerFile << "Maximum axial difference mm: " << m_maxAxialDiffmm << endl;
1121  // Maximum number of lines per event (write it only if higher than 1
1122  if (m_maxNumberOfLinesPerEvent >1) headerFile << "Maximum number of lines per event: " << m_maxNumberOfLinesPerEvent << endl;
1123 
1124  // Calibration factor
1125  headerFile << "Calibration factor: " << m_calibrationFactor << endl;
1126  // Isotope
1127  headerFile << "Isotope: " << m_isotope << endl;
1128  // TOF info
1129  if (m_TOFInfoFlag)
1130  {
1131  // Flag
1132  headerFile << "TOF information flag: 1" << endl;
1133  // TOF resolution
1134  headerFile << "TOF resolution (ps): " << m_resolutionTOF << endl;
1135  // List-mode TOF measurement range
1136  if (m_dataMode == MODE_LIST) headerFile << "List TOF measurement range (ps): " << m_TOFMeasurementRange << endl;
1137  // Histogram
1138  else if (m_dataMode == MODE_HISTOGRAM)
1139  {
1140  // Number of TOF bins
1141  headerFile << "Histo TOF number of bins: " << m_nbBinsTOF << endl;
1142  // TOF bin size
1143  headerFile << "Histo TOF bin size (ps): " << m_binSizeTOF << endl;
1144  }
1145  }
1146  // Event kind
1147  if (m_eventKindFlag)
1148  {
1149  // Error if histogram data are used
1151  {
1152  Cerr("***** iDataFilePET::WriteHeader -> Request for writing event type in histogram mode (this field is specific to list-mode) !" << endl);
1153  return 1;
1154  }
1155  else headerFile << "Coincidence kind flag: " << m_eventKindFlag << endl;
1156  }
1157  // Correction flags
1158  if (m_atnCorrectionFlag)
1159  headerFile << "Attenuation correction flag: " << m_atnCorrectionFlag << endl;
1161  headerFile << "Normalization correction flag: " << m_normCorrectionFlag << endl;
1163  headerFile << "Scatter correction flag: " << m_scatCorrectionFlag << endl;
1165  headerFile << "Random correction flag: " << m_randCorrectionFlag << endl;
1166  // Close file
1167  headerFile.close();
1168  // End
1169  return 0;
1170 }
1171 
1172 // =====================================================================
1173 // ---------------------------------------------------------------------
1174 // ---------------------------------------------------------------------
1175 // =====================================================================
1176 
1178 {
1179  // Verbose
1181  if (m_verbose>=VERBOSE_DETAIL) Cout("PROJ_GetScannerSpecificParameters() ..."<< endl);
1182  // Get pointers to SPECT specific parameters in scanner
1183  if (sScannerManager::GetInstance()->PROJ_GetPETSpecificParameters(&m_maxAxialDiffmm) )
1184  {
1185  Cerr("***** iDataFilePET::PROJ_GetScannerSpecificParameters() -> An error occurred while trying to get PET geometric parameters from the scanner object !" << endl);
1186  return 1;
1187  }
1188  // End
1189  return 0;
1190 }
1191 
1192 // =====================================================================
1193 // ---------------------------------------------------------------------
1194 // ---------------------------------------------------------------------
1195 // =====================================================================
int WriteEvent(vEvent *ap_Event, int a_th)
Write event according to the chosen type of data.
bool m_randCorrectionFlag
This class is designed to be a mother virtual class for DataFile.
Definition: vDataFile.hh:103
bool GetNormCorrectionFlag()
Simply return m_normCorrectionFlag.
int SetPETIsotope(int a_bed, const string &a_isotope)
Set the PET isotope for the provided bed.
static sScannerManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
#define VERBOSE_DEBUG_EVENT
bool GetTOFInfoFlag()
#define MODE_HISTOGRAM
Definition: vDataFile.hh:59
#define TYPE_PET
Definition: vDataFile.hh:74
void SetNbLines(uint16_t a_value)
Set the number of lines of the Event.
Definition: vEvent.hh:154
int64_t m_sizeEvent
Definition: vDataFile.hh:598
bool mp_POIDirectionFlag[3]
Definition: vDataFile.hh:594
string m_dataFileName
Definition: vDataFile.hh:578
FLTNB GetFrameTimeStartInSec(int a_bed, int a_frame)
Get the frame time start for the given bed, in seconds as a FLTNB.
#define MODE_LIST
Definition: vDataFile.hh:57
int PROJ_GetScannerSpecificParameters()
Get PET specific parameters for projections from the scanner object, through the scannerManager.
#define FLTNB
Definition: gVariables.hh:81
#define MODE_NORMALIZATION
Definition: vDataFile.hh:61
FLTNB m_binSizeTOF
bool m_ignoreNormCorrectionFlag
bool m_ignoreTOFFlag
void DescribeSpecific()
Implementation of the pure virtual eponym function that simply prints info about the datafile...
FLTNB GetEventScatRate(int a_bin)
FLTNB m_maxAxialDiffmm
uint32_t GetID2(int a_line)
Definition: vEvent.hh:108
#define VERBOSE_DETAIL
int CheckSpecificConsistencyWithAnotherDataFile(vDataFile *ap_DataFile)
Check consistency between 'this' and the provided datafile, for specific characteristics.
FLTNB m_durationInSec
Definition: vDataFile.hh:584
iDataFilePET()
iDataFilePET constructor. Initialize the member variables to their default values.
Definition: iDataFilePET.cc:39
bool GetIgnoreAttnCorrectionFlag()
Get the boolean that says if the attenuation correction is ignored or not.
int ReadSpecificInfoInHeader(bool a_affectQuantificationFlag)
Read through the header file and gather specific PET information.
Definition: iDataFilePET.cc:76
bool IsInitialized()
Returns true if the object has been initialized.
string m_headerFileName
Definition: vDataFile.hh:577
int CheckSpecificParameters()
Check parameters specific to PET data.
int m_dataMode
Definition: vDataFile.hh:580
static sOutputManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
vEvent * GetEventSpecific(char *ap_buffer, int a_th)
Read an event from the position pointed by 'ap_buffer', parse the generic or modality-specific inform...
bool m_ignoreScatCorrectionFlag
int ComputeSizeEvent()
Computation of the size of each event according to the mandatory/optional correction fields...
Declaration of class iDataFilePET.
int m_dataSpec
Definition: vDataFile.hh:582
FLTNB m_startTimeInSec
Definition: vDataFile.hh:583
bool m_eventKindFlag
string GetFileFromPath(const string &a_pathToFile)
Simply return the file from a path string passed in parameter.
Definition: gOptions.cc:1141
string GetIsotope()
bool GetIgnoreNormCorrectionFlag()
Get the boolean that says if the normalization correction is ignored or not.
bool GetIgnoreRandCorrectionFlag()
Get the boolean that says if the random correction is ignored or not.
#define FLTNBDATA
Definition: gVariables.hh:87
FLTNB GetTOFMeasurement()
void SetAttenuationCorrectionFactor(FLTNBDATA a_value)
Cast the FLTNBDATA value passed as parameter in FLTNB, and set it to the attenuation correction facto...
Definition: iEventNorm.hh:89
Inherit from iEventPET. Class for PET list-mode events.
FLTNB GetEventValue(int a_bin)
#define Cerr(MESSAGE)
const string & GetPathName()
Inherit from iEventPET. Class for PET histogram mode events.
#define VERBOSE_DEBUG_LIGHT
bool m_ignorePOIFlag
Definition: vDataFile.hh:593
~iDataFilePET()
iDataFilePET destructor.
Definition: iDataFilePET.cc:69
FLTNB GetTOFResolution()
Singleton class that manages output writing on disk (images, sinograms, etc). It also manages loggi...
FLTNB m_TOFMeasurementRange
FLTNB m_resolutionTOF
int ReadDataASCIIFile(const string &a_file, const string &a_keyword, T *ap_return, int a_nbElts, bool a_mandatoryFlag)
Look for "a_nbElts" elts in the "a_file" file matching the "a_keyword" string passed as parameter a...
Definition: gOptions.cc:123
bool m_ignoreAttnCorrectionFlag
int m_dataType
Definition: vDataFile.hh:581
int64_t m_nbEvents
Definition: vDataFile.hh:579
bool GetScatCorrectionFlag()
Simply return m_scatCorrectionFlag.
int SetSpecificParametersFrom(vDataFile *ap_DataFile)
Initialize all parameters specific to PET from the provided datafile.
oImageDimensionsAndQuantification * mp_ID
Definition: vDataFile.hh:573
uint16_t GetMaxNumberOfLinesPerEvent()
fstream ** m2p_dataFile
Definition: vDataFile.hh:599
bool GetEventKindFlag()
Simply return m_eventKindFlag.
#define VERBOSE_NORMAL
uint8_t GetKind()
FLTNB GetNormFactor()
Definition: iEventPET.hh:85
int WriteHistoEvent(iEventHistoPET *ap_Event, int a_th)
Write a PET histogram event.
int WriteHeader()
Generate a header file according to the data output information.
int m_bedIndex
Definition: vDataFile.hh:586
int PROJ_InitFile()
Initialize the fstream objets for output writing as well as some other variables specific to the Proj...
#define KEYWORD_OPTIONAL
Definition: gOptions.hh:50
FLTNB * GetPOI2()
Inherit from vEvent. Used for normalization events for sensitivity computation.
Definition: iEventNorm.hh:43
int PrepareDataFile()
Store different kind of information inside arrays (data relative to specific correction as well as ba...
FLTNB mp_POIResolution[3]
Definition: vDataFile.hh:595
const string & GetBaseName()
uint16_t GetNbLines()
Definition: vEvent.hh:94
void SetTimeInMs(uint32_t a_value)
Set the timestamp of the Event.
Definition: vEvent.hh:139
FLTNB GetMaxAxialDiffmm()
Mother class for the Event objects.
Definition: vEvent.hh:43
FLTNB GetTOFBinSize()
int GetDataMode()
Definition: vDataFile.hh:300
string m_isotope
int GetNbTimeFrames()
Get the number of time frames.
FLTNB GetAtnCorrFactor()
Definition: iEventPET.hh:91
FLTNB GetTOFMeasurementRange()
vEvent ** m2p_BufferEvent
Definition: vDataFile.hh:600
bool GetRandCorrectionFlag()
Simply return m_randCorrectionFlag.
bool m_POIInfoFlag
Definition: vDataFile.hh:592
#define SPEC_EMISSION
Definition: vDataFile.hh:91
uint32_t GetID1(int a_line)
Definition: vEvent.hh:101
FLTNB * GetPOI1()
#define SPEED_OF_LIGHT
Definition: gVariables.hh:106
int GetNbThreadsForProjection()
Get the number of threads used for projections.
bool GetIgnoreScatCorrectionFlag()
Get the boolean that says if the scatter correction is ignored or not.
int GetNbTOFBins()
#define DEBUG_VERBOSE(IGNORED1, IGNORED2)
int CheckFileSizeConsistency()
This function is implemented in child classes Check if file size is consistent. ...
#define Cout(MESSAGE)
FLTNB m_calibrationFactor
Definition: vDataFile.hh:585
FLTNB GetEventRdmRate()
Definition: iEventPET.hh:79
FLTNB GetFrameDurationInSec(int a_bed, int a_frame)
Get the frame duration for the given bed, in seconds as a FLTNB.
bool m_atnCorrectionFlag
bool m_normCorrectionFlag
uint16_t m_maxNumberOfLinesPerEvent
bool GetAtnCorrectionFlag()
Simply return m_atnCorrectionFlag.
Inherit from vDataFile. Class that manages the reading of a PET input file (header + data)...
Definition: iDataFilePET.hh:44
bool m_scatCorrectionFlag
FLTNB GetEventScatRate(int a_bin)
int m_verbose
Definition: vDataFile.hh:574
uint32_t GetTimeInMs()
Definition: vEvent.hh:88
bool m_ignoreRandCorrectionFlag
int WriteListEvent(iEventListPET *ap_Event, int a_th=0)
Write a PET list-mode event.