/
itkScancoImageIO.h
374 lines (312 loc) · 9.78 KB
/
itkScancoImageIO.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* http://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
/*=========================================================================
Program: DICOM for VTK
Copyright (c) 2015 David Gobbi
All rights reserved.
See Copyright.txt or http://dgobbi.github.io/bsd3.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
/**
* David's notes from vtkScancoCTReader
* Read SCANCO ISQ and AIM medical image files
*
* This class reads ISQ and AIM files, which are used for high-resolution
* computed tomography. The information that it provides uses different
* units as compared to the original files: all distances are given in
* millimeters (instead of micrometers), times are given in milliseconds
* (instead of microseconds), voltage and current given in kV and mA
* (instead of volts and microamps). If the scanner was calibrated, then
* the data values can be converted to calibrated units. To convert
* to linear attenuation coefficients [cm^-1], simply divide the data
* values by the MuScaling. To convert to density values, multiply
* the data values by the m_RescaleSlope and add the m_RescaleIntercept.
* To convert to Hounsfield units, multiply by 1000/(MuScaling*m_MuWater)
* and subtract 1000.
*
* Created at the Calgary Image Processing and Analysis Centre (CIPAC).
*/
#ifndef itkScancoImageIO_h
#define itkScancoImageIO_h
#include "IOScancoExport.h"
#include <fstream>
#include "itkImageIOBase.h"
namespace itk
{
/** \class ScancoImageIO
*
* \brief Read Scanco image file formats.
*
* Many methods are based off vtkScancoCTReader in vtk-dicom by David Gobbi
*
* \ingroup IOFilters
* \ingroup IOScanco
*/
class IOScanco_EXPORT ScancoImageIO : public ImageIOBase
{
public:
ITK_DISALLOW_COPY_AND_MOVE(ScancoImageIO);
/** Standard class typedefs. */
using Self = ScancoImageIO;
using Superclass = ImageIOBase;
using Pointer = SmartPointer<Self>;
using ConstPointer = SmartPointer<const Self>;
/** Method for creation through the object factory. */
itkNewMacro(Self);
/** Run-time type information (and related methods). */
itkTypeMacro(ScancoImageIO, ImageIOBase);
/** The different types of ImageIO's can support data of varying
* dimensionality. For example, some file formats are strictly 2D
* while others can support 2D, 3D, or even n-D. This method returns
* true/false as to whether the ImageIO can support the dimension
* indicated. */
bool
SupportsDimension(unsigned long dimension) override
{
if (dimension == 3)
{
return true;
}
return false;
}
/*-------- This part of the interfaces deals with reading data. ----- */
/** Determine the file type. Returns true if this ImageIO can read the
* file specified. */
bool
CanReadFile(const char *) override;
/** Set the spacing and dimension information for the set filename. */
void
ReadImageInformation() override;
/** Reads the data from disk into the memory buffer provided. */
void
Read(void * buffer) override;
/*-------- This part of the interfaces deals with writing data. ----- */
/** Determine the file type. Returns true if this ImageIO can write the
* file specified. */
bool
CanWriteFile(const char *) override;
/** Set the spacing and dimension information for the set filename. */
void
WriteImageInformation() override;
/** Writes the data to disk from the memory buffer provided. Make sure
* that the IORegions has been set properly. */
void
Write(const void * buffer) override;
bool
CanStreamRead() override
{
return false;
}
bool
CanStreamWrite() override
{
return false;
}
/** Get a string that states the version of the file header.
* Max size: 16 characters. */
const char *
GetVersion() const
{
return this->m_Version;
}
void
SetVersion(const char * version)
{
strncpy(this->m_Version, version, 18);
this->Modified();
}
itkGetConstMacro(PatientIndex, int);
itkSetMacro(PatientIndex, int);
itkGetConstMacro(ScannerID, int);
itkSetMacro(ScannerID, int);
itkGetConstMacro(SliceThickness, double);
itkSetMacro(SliceThickness, double);
itkGetConstMacro(SliceIncrement, double);
itkSetMacro(SliceIncrement, double);
itkGetConstMacro(StartPosition, double);
itkSetMacro(StartPosition, double);
/** Set / Get the minimum and maximum values */
const double *
GetDataRange() const
{
return this->m_DataRange;
}
void
SetDataRange(const double * dataRange)
{
this->m_DataRange[0] = dataRange[0];
this->m_DataRange[1] = dataRange[1];
}
itkGetConstMacro(MuScaling, double);
itkSetMacro(MuScaling, double);
itkGetConstMacro(NumberOfSamples, int);
itkSetMacro(NumberOfSamples, int);
itkGetConstMacro(NumberOfProjections, int);
itkSetMacro(NumberOfProjections, int);
itkGetConstMacro(ScanDistance, double);
itkSetMacro(ScanDistance, double);
itkGetConstMacro(ScannerType, int);
itkSetMacro(ScannerType, int);
itkGetConstMacro(SampleTime, double);
itkSetMacro(SampleTime, double);
itkGetConstMacro(MeasurementIndex, int);
itkSetMacro(MeasurementIndex, int);
itkGetConstMacro(Site, int);
itkSetMacro(Site, int);
itkGetConstMacro(ReferenceLine, int);
itkSetMacro(ReferenceLine, int);
itkGetConstMacro(ReconstructionAlg, int);
itkSetMacro(ReconstructionAlg, int);
/** Get a string that states patient name.
* Max size: 40 characters. */
const char *
GetPatientName() const
{
return this->m_PatientName;
}
void
SetPatientName(const char * patientName)
{
strncpy(this->m_PatientName, patientName, 42);
this->Modified();
}
const char *
GetCreationDate() const
{
return this->m_CreationDate;
}
void
SetCreationDate(const char * creationDate)
{
strncpy(this->m_CreationDate, creationDate, 32);
this->Modified();
}
const char *
GetModificationDate() const
{
return this->m_ModificationDate;
}
void
SetModificationDate(const char * modificationDate)
{
strncpy(this->m_ModificationDate, modificationDate, 32);
this->Modified();
}
itkGetConstMacro(Energy, double);
itkSetMacro(Energy, double);
itkGetConstMacro(Intensity, double);
itkSetMacro(Intensity, double);
protected:
ScancoImageIO();
~ScancoImageIO() override;
void
PrintSelf(std::ostream & os, Indent indent) const override;
private:
/** Check the file header to see what type of file it is.
*
* Return values are: 0 if unrecognized, 1 if ISQ/RAD,
* 2 if AIM 020, 3 if AIM 030.
*/
int
CheckVersion(const char header[16]);
/** Convert char data to 32-bit int (little-endian). */
static int
DecodeInt(const void * data);
/** Convert 32-bit int (little-endian) to char data. */
static void
EncodeInt(int data, void * target);
/** Convert char data to float (single precision). */
static float
DecodeFloat(const void * data);
/** Convert char data to float (double precision). */
static double
DecodeDouble(const void * data);
//! Convert a VMS timestamp to a calendar date.
void
DecodeDate(const void * data,
int & year,
int & month,
int & day,
int & hour,
int & minute,
int & second,
int & millis);
//! Convert the current calendar date to a VMS timestamp and store in target
void
EncodeDate(void * target);
//! Strip a string by removing trailing whitespace.
/*!
* The dest must have a size of at least l+1.
*/
static void
StripString(char * dest, const char * source, size_t length);
static void
PadString(char * dest, const char * source, size_t length);
void
InitializeHeader();
int
ReadISQHeader(std::ifstream * file, unsigned long bytesRead);
int
ReadAIMHeader(std::ifstream * file, unsigned long bytesRead);
void
PopulateMetaDataDictionary();
void
WriteISQHeader(std::ofstream * file);
// Header information
char m_Version[18];
char m_PatientName[42];
int m_PatientIndex;
int m_ScannerID;
char m_CreationDate[32];
char m_ModificationDate[32];
int ScanDimensionsPixels[3];
double ScanDimensionsPhysical[3];
double m_SliceThickness;
double m_SliceIncrement;
double m_StartPosition;
double m_EndPosition;
double m_ZPosition;
double m_DataRange[2];
double m_MuScaling;
int m_NumberOfSamples;
int m_NumberOfProjections;
double m_ScanDistance;
double m_SampleTime;
int m_ScannerType;
int m_MeasurementIndex;
int m_Site;
int m_ReconstructionAlg;
double m_ReferenceLine;
double m_Energy;
double m_Intensity;
int m_RescaleType;
char m_RescaleUnits[18];
char m_CalibrationData[66];
double m_RescaleSlope;
double m_RescaleIntercept;
double m_MuWater;
char * m_RawHeader;
// The compression mode, if any.
int m_Compression;
bool m_HeaderInitialized = false;
SizeValueType m_HeaderSize{ 0 };
};
} // end namespace itk
#endif // itkScancoImageIO_h