1 /*****************************************************************************
4 ** Name: projections.cpp Projection data classes
5 ** Programmer: Kevin Rosenberg
6 ** Date Started: Aug 84
8 ** This is part of the CTSim program
9 ** Copyright (c) 1983-2001 Kevin Rosenberg
11 ** $Id: projections.cpp,v 1.73 2001/03/30 19:17:32 kevin Exp $
13 ** This program is free software; you can redistribute it and/or modify
14 ** it under the terms of the GNU General Public License (version 2) as
15 ** published by the Free Software Foundation.
17 ** This program is distributed in the hope that it will be useful,
18 ** but WITHOUT ANY WARRANTY; without even the implied warranty of
19 ** MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20 ** GNU General Public License for more details.
22 ** You should have received a copy of the GNU General Public License
23 ** along with this program; if not, write to the Free Software
24 ** Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25 ******************************************************************************/
29 const kuint16 Projections::m_signature = ('P'*256 + 'J');
31 const int Projections::POLAR_INTERP_INVALID = -1;
32 const int Projections::POLAR_INTERP_NEAREST = 0;
33 const int Projections::POLAR_INTERP_BILINEAR = 1;
34 const int Projections::POLAR_INTERP_BICUBIC = 2;
36 const char* const Projections::s_aszInterpName[] =
43 const char* const Projections::s_aszInterpTitle[] =
50 const int Projections::s_iInterpCount = sizeof(s_aszInterpName) / sizeof(char*);
55 * Projections Constructor for projections matrix storage
58 * proj = projections_create (filename, nView, nDet)
59 * Projections& proj Allocated projections structure & matrix
60 * int nView Number of rotated view
61 * int nDet Number of detectors
65 Projections::Projections (const Scanner& scanner)
68 initFromScanner (scanner);
72 Projections::Projections (const int nView, const int nDet)
78 Projections::Projections (void)
84 Projections::~Projections (void)
90 Projections::convertInterpNameToID (const char* const interpName)
92 int interpID = POLAR_INTERP_INVALID;
94 for (int i = 0; i < s_iInterpCount; i++)
95 if (strcasecmp (interpName, s_aszInterpName[i]) == 0) {
104 Projections::convertInterpIDToName (const int interpID)
106 static const char *interpName = "";
108 if (interpID >= 0 && interpID < s_iInterpCount)
109 return (s_aszInterpName[interpID]);
115 Projections::convertInterpIDToTitle (const int interpID)
117 static const char *interpTitle = "";
119 if (interpID >= 0 && interpID < s_iInterpCount)
120 return (s_aszInterpTitle[interpID]);
122 return (interpTitle);
128 Projections::init (const int nView, const int nDet)
130 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
135 time_t t = time (NULL);
136 tm* lt = localtime (&t);
137 m_year = lt->tm_year;
138 m_month = lt->tm_mon;
140 m_hour = lt->tm_hour;
141 m_minute = lt->tm_min;
142 m_second = lt->tm_sec;
146 Projections::initFromScanner (const Scanner& scanner)
148 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
150 init (scanner.nView(), scanner.nDet());
152 m_rotInc = scanner.rotInc();
153 m_detInc = scanner.detInc();
154 m_detStart = scanner.detStart();
155 m_geometry = scanner.geometry();
156 m_dFocalLength = scanner.focalLength();
157 m_dSourceDetectorLength = scanner.sourceDetectorLength();
158 m_dViewDiameter = scanner.viewDiameter();
160 m_dFanBeamAngle = scanner.fanBeamAngle();
164 Projections::setNView (int nView) // used by MPI to reduce # of views
167 init (nView, m_nDet);
174 Projections::newProjData (void)
177 sys_error(ERR_WARNING, "m_projData != NULL [newProjData]");
179 if (m_nView > 0 && m_nDet) {
180 m_projData = new DetectorArray* [m_nView];
182 for (int i = 0; i < m_nView; i++)
183 m_projData[i] = new DetectorArray (m_nDet);
189 * projections_free Free memory allocated to projections
192 * projections_free(proj)
193 * Projections& proj Projectionss to be deallocated
197 Projections::deleteProjData (void)
199 if (m_projData != NULL) {
200 for (int i = 0; i < m_nView; i++)
201 delete m_projData[i];
210 * Projections::headerWwrite Write data header for projections file
215 Projections::headerWrite (fnetorderstream& fs)
217 kuint16 _hsize = m_headerSize;
218 kuint16 _signature = m_signature;
219 kuint32 _nView = m_nView;
220 kuint32 _nDet = m_nDet;
221 kuint32 _geom = m_geometry;
222 kuint16 _remarksize = m_remark.length();
223 kuint16 _year = m_year;
224 kuint16 _month = m_month;
225 kuint16 _day = m_day;
226 kuint16 _hour = m_hour;
227 kuint16 _minute = m_minute;
228 kuint16 _second = m_second;
230 kfloat64 _calcTime = m_calcTime;
231 kfloat64 _rotStart = m_rotStart;
232 kfloat64 _rotInc = m_rotInc;
233 kfloat64 _detStart = m_detStart;
234 kfloat64 _detInc = m_detInc;
235 kfloat64 _viewDiameter = m_dViewDiameter;
236 kfloat64 _focalLength = m_dFocalLength;
237 kfloat64 _sourceDetectorLength = m_dSourceDetectorLength;
238 kfloat64 _fanBeamAngle = m_dFanBeamAngle;
244 fs.writeInt16 (_hsize);
245 fs.writeInt16 (_signature);
246 fs.writeInt32 (_nView);
247 fs.writeInt32 (_nDet);
248 fs.writeInt32 (_geom);
249 fs.writeFloat64 (_calcTime);
250 fs.writeFloat64 (_rotStart);
251 fs.writeFloat64 (_rotInc);
252 fs.writeFloat64 (_detStart);
253 fs.writeFloat64 (_detInc);
254 fs.writeFloat64 (_viewDiameter);
255 fs.writeFloat64 (_focalLength);
256 fs.writeFloat64 (_sourceDetectorLength);
257 fs.writeFloat64 (_fanBeamAngle);
258 fs.writeInt16 (_year);
259 fs.writeInt16 (_month);
260 fs.writeInt16 (_day);
261 fs.writeInt16 (_hour);
262 fs.writeInt16 (_minute);
263 fs.writeInt16 (_second);
264 fs.writeInt16 (_remarksize);
265 fs.write (m_remark.c_str(), _remarksize);
267 m_headerSize = fs.tellp();
268 _hsize = m_headerSize;
270 fs.writeInt16 (_hsize);
278 * projections_read_header Read data header for projections file
282 Projections::headerRead (fnetorderstream& fs)
284 kuint16 _hsize, _signature, _year, _month, _day, _hour, _minute, _second, _remarksize = 0;
285 kuint32 _nView, _nDet, _geom;
286 kfloat64 _calcTime, _rotStart, _rotInc, _detStart, _detInc, _focalLength, _sourceDetectorLength, _viewDiameter, _fanBeamAngle;
292 fs.readInt16 (_hsize);
293 fs.readInt16 (_signature);
294 fs.readInt32 (_nView);
295 fs.readInt32 (_nDet);
296 fs.readInt32 (_geom);
297 fs.readFloat64 (_calcTime);
298 fs.readFloat64 (_rotStart);
299 fs.readFloat64 (_rotInc);
300 fs.readFloat64 (_detStart);
301 fs.readFloat64 (_detInc);
302 fs.readFloat64 (_viewDiameter);
303 fs.readFloat64 (_focalLength);
304 fs.readFloat64 (_sourceDetectorLength);
305 fs.readFloat64 (_fanBeamAngle);
306 fs.readInt16 (_year);
307 fs.readInt16 (_month);
309 fs.readInt16 (_hour);
310 fs.readInt16 (_minute);
311 fs.readInt16 (_second);
312 fs.readInt16 (_remarksize);
315 sys_error (ERR_SEVERE, "Error reading header information , _remarksize=%d [projections_read_header]", _remarksize);
319 if (_signature != m_signature) {
320 sys_error (ERR_SEVERE, "File %s does not have a valid projection file signature", m_filename.c_str());
324 char* pszRemarkStorage = new char [_remarksize+1];
325 fs.read (pszRemarkStorage, _remarksize);
327 sys_error (ERR_SEVERE, "Error reading remark, _remarksize = %d", _remarksize);
330 pszRemarkStorage[_remarksize] = 0;
331 m_remark = pszRemarkStorage;
332 delete pszRemarkStorage;
334 off_t _hsizeread = fs.tellg();
335 if (!fs || _hsizeread != _hsize) {
336 sys_error (ERR_WARNING, "File header size read %ld != file header size stored %ld [read_projections_header]\n_remarksize=%ld", (long int) _hsizeread, _hsize, _remarksize);
340 m_headerSize = _hsize;
344 m_calcTime = _calcTime;
345 m_rotStart = _rotStart;
347 m_detStart = _detStart;
349 m_dFocalLength = _focalLength;
350 m_dSourceDetectorLength = _sourceDetectorLength;
351 m_dViewDiameter = _viewDiameter;
352 m_dFanBeamAngle = _fanBeamAngle;
360 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
361 m_label.setLabelString (m_remark);
362 m_label.setCalcTime (m_calcTime);
363 m_label.setDateTime (m_year, m_month, m_day, m_hour, m_minute, m_second);
369 Projections::read (const std::string& filename)
371 return read (filename.c_str());
376 Projections::read (const char* filename)
378 m_filename = filename;
380 frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary);
382 frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary | std::ios::nocreate);
388 if (! headerRead (fileRead))
394 for (int i = 0; i < m_nView; i++) {
395 if (! detarrayRead (fileRead, *m_projData[i], i))
405 Projections::copyViewData (const std::string& filename, std::ostream& os, int startView, int endView)
407 return copyViewData (filename.c_str(), os, startView, endView);
411 Projections::copyViewData (const char* const filename, std::ostream& os, int startView, int endView)
413 frnetorderstream is (filename, std::ios::in | std::ios::binary);
414 kuint16 sizeHeader, signature;
415 kuint32 _nView, _nDet;
419 sys_error (ERR_SEVERE, "Unable to read projection file %s", filename);
423 is.readInt16 (sizeHeader);
424 is.readInt16 (signature);
425 is.readInt32 (_nView);
426 is.readInt32 (_nDet);
430 if (signature != m_signature) {
431 sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
437 if (startView > nView - 1)
439 if (endView < 0 || endView > nView - 1)
442 if (startView > endView) { // swap if start > end
443 int tempView = endView;
445 startView = tempView;
448 int sizeView = 8 /* view_angle */ + 4 /* nDet */ + (4 * nDet);
449 unsigned char* pViewData = new unsigned char [sizeView];
451 for (int i = startView; i <= endView; i++) {
452 is.seekg (sizeHeader + i * sizeView);
453 is.read (reinterpret_cast<char*>(pViewData), sizeView);
454 os.write (reinterpret_cast<char*>(pViewData), sizeView);
455 if (is.fail() || os.fail())
461 sys_error (ERR_SEVERE, "Error reading projection file");
463 sys_error (ERR_SEVERE, "Error writing projection file");
465 return (! (is.fail() | os.fail()));
469 Projections::copyHeader (const std::string& filename, std::ostream& os)
471 return copyHeader (filename.c_str(), os);
475 Projections::copyHeader (const char* const filename, std::ostream& os)
477 frnetorderstream is (filename, std::ios::in | std::ios::binary);
478 kuint16 sizeHeader, signature;
479 is.readInt16 (sizeHeader);
480 is.readInt16 (signature);
482 if (signature != m_signature) {
483 sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
487 unsigned char* pHdrData = new unsigned char [sizeHeader];
488 is.read (reinterpret_cast<char*>(pHdrData), sizeHeader);
490 sys_error (ERR_SEVERE, "Error reading header");
494 os.write (reinterpret_cast<char*>(pHdrData), sizeHeader);
496 sys_error (ERR_SEVERE, "Error writing header");
504 Projections::write (const std::string& filename)
506 return write (filename.c_str());
510 Projections::write (const char* filename)
512 frnetorderstream fs (filename, std::ios::out | std::ios::binary | std::ios::trunc | std::ios::ate);
513 m_filename = filename;
515 sys_error (ERR_SEVERE, "Error opening file %s for output [projections_create]", filename);
519 if (! headerWrite (fs))
522 if (m_projData != NULL) {
523 for (int i = 0; i < m_nView; i++) {
524 if (! detarrayWrite (fs, *m_projData[i], i))
537 * detarrayRead Read a Detector Array structure from the disk
540 * detarrayRead (proj, darray, view_num)
541 * DETARRAY *darray Detector array storage location to be filled
542 * int view_num View number to read
546 Projections::detarrayRead (fnetorderstream& fs, DetectorArray& darray, const int iview)
548 const int detval_bytes = darray.nDet() * sizeof(kfloat32);
549 const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
550 const int view_bytes = detheader_bytes + detval_bytes;
551 const off_t start_data = m_headerSize + (iview * view_bytes);
552 DetectorValue* detval_ptr = darray.detValues();
556 fs.seekg (start_data);
558 fs.readFloat64 (view_angle);
560 darray.setViewAngle (view_angle);
561 // darray.setNDet ( nDet);
563 for (unsigned int i = 0; i < nDet; i++) {
565 fs.readFloat32 (detval);
566 detval_ptr[i] = detval;
576 * detarrayWrite Write detector array data to the disk
579 * detarrayWrite (darray, view_num)
580 * DETARRAY *darray Detector array data to be written
581 * int view_num View number to write
584 * This routine writes the detarray data from the disk sequentially to
585 * the file that was opened with open_projections(). Data is written in
590 Projections::detarrayWrite (fnetorderstream& fs, const DetectorArray& darray, const int iview)
592 const int detval_bytes = darray.nDet() * sizeof(float);
593 const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
594 const int view_bytes = detheader_bytes + detval_bytes;
595 const off_t start_data = m_headerSize + (iview * view_bytes);
596 const DetectorValue* const detval_ptr = darray.detValues();
597 kfloat64 view_angle = darray.viewAngle();
598 kuint32 nDet = darray.nDet();
600 fs.seekp (start_data);
602 sys_error (ERR_SEVERE, "Error seeking detectory array [detarrayWrite]");
606 fs.writeFloat64 (view_angle);
607 fs.writeInt32 (nDet);
609 for (unsigned int i = 0; i < nDet; i++) {
610 kfloat32 detval = detval_ptr[i];
611 fs.writeFloat32 (detval);
621 * printProjectionData Print projections data
624 * printProjectionData ()
628 Projections::printProjectionData ()
630 printProjectionData (0, nView() - 1);
634 Projections::printProjectionData (int startView, int endView)
636 printf("Projections Data\n\n");
637 printf("Description: %s\n", m_remark.c_str());
638 printf("Geometry: %s\n", Scanner::convertGeometryIDToName (m_geometry));
639 printf("nView = %8d nDet = %8d\n", m_nView, m_nDet);
640 printf("focalLength = %8.4f ViewDiameter = %8.4f\n", m_dFocalLength, m_dViewDiameter);
641 printf("fanBeamAngle= %8.4f SourceDetector = %8.4f\n", convertRadiansToDegrees(m_dFanBeamAngle), m_dSourceDetectorLength);
642 printf("rotStart = %8.4f rotInc = %8.4f\n", m_rotStart, m_rotInc);
643 printf("detStart = %8.4f detInc = %8.4f\n", m_detStart, m_detInc);
644 if (m_projData != NULL) {
648 endView = m_nView - 1;
649 if (startView > m_nView - 1)
650 startView = m_nView - 1;
651 if (endView > m_nView - 1)
652 endView = m_nView - 1;
653 for (int ir = startView; ir <= endView - 1; ir++) {
654 printf("View %d: angle %f\n", ir, m_projData[ir]->viewAngle());
655 DetectorValue* detval = m_projData[ir]->detValues();
656 for (int id = 0; id < m_projData[ir]->nDet(); id++)
657 printf("%8.4f ", detval[id]);
664 Projections::printScanInfo (std::ostringstream& os) const
666 os << "Number of detectors: " << m_nDet << "\n";
667 os << "Number of views: " << m_nView<< "\n";
668 os << "Description: " << m_remark.c_str()<< "\n";
669 os << "Geometry: " << Scanner::convertGeometryIDToName (m_geometry)<< "\n";
670 os << "Focal Length: " << m_dFocalLength<< "\n";
671 os << "Source Detector Length: " << m_dSourceDetectorLength << "\n";
672 os << "View Diameter: " << m_dViewDiameter<< "\n";
673 os << "Fan Beam Angle: " << convertRadiansToDegrees(m_dFanBeamAngle) << "\n";
674 os << "detStart: " << m_detStart<< "\n";
675 os << "detInc: " << m_detInc<< "\n";
676 os << "rotStart: " << m_rotStart<< "\n";
677 os << "rotInc: " << m_rotInc<< "\n";
682 Projections::convertPolar (ImageFile& rIF, int iInterpolationID)
684 unsigned int nx = rIF.nx();
685 unsigned int ny = rIF.ny();
686 ImageFileArray v = rIF.getArray();
687 ImageFileArray vImag = rIF.getImaginaryArray();
689 if (! v || nx == 0 || ny == 0)
692 Projections* pProj = this;
693 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
694 pProj = interpolateToParallel();
696 Array2d<double> adView (nx, ny);
697 Array2d<double> adDet (nx, ny);
698 double** ppdView = adView.getArray();
699 double** ppdDet = adDet.getArray();
701 std::complex<double>** ppcDetValue = new std::complex<double>* [pProj->m_nView];
703 for (iView = 0; iView < pProj->m_nView; iView++) {
704 ppcDetValue[iView] = new std::complex<double> [pProj->m_nDet];
705 DetectorValue* detval = pProj->getDetectorArray (iView).detValues();
706 for (unsigned int iDet = 0; iDet < pProj->m_nDet; iDet++)
707 ppcDetValue[iView][iDet] = std::complex<double>(detval[iDet], 0);
710 pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet, pProj->m_nDet, 1., pProj->m_detInc);
712 pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet,
713 pProj->m_nDet, iInterpolationID);
715 for (iView = 0; iView < pProj->m_nView; iView++)
716 delete [] ppcDetValue[iView];
717 delete [] ppcDetValue;
719 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
727 Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad)
730 rIF.arrayDataClear();
733 unsigned int nx = rIF.nx();
734 unsigned int ny = rIF.ny();
735 ImageFileArray v = rIF.getArray();
736 if (! rIF.isComplex())
737 rIF.convertRealToComplex();
738 ImageFileArray vImag = rIF.getImaginaryArray();
740 if (! v || nx == 0 || ny == 0)
743 Projections* pProj = this;
744 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
745 pProj = interpolateToParallel();
748 // int iInterpDet = pProj->m_nDet;
749 int iNumInterpDetWithZeros = ProcessSignal::addZeropadFactor (iInterpDet, iZeropad);
751 double dZeropadRatio = static_cast<double>(iNumInterpDetWithZeros) / static_cast<double>(iInterpDet);
753 fftw_plan plan = fftw_create_plan (iNumInterpDetWithZeros, FFTW_FORWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
755 fftw_complex* pcIn = new fftw_complex [iNumInterpDetWithZeros];
756 std::complex<double>** ppcDetValue = new std::complex<double>* [pProj->m_nView];
757 double dInterpScale = (pProj->m_nDet-1) / static_cast<double>(iInterpDet-1) / SQRT2;
759 double dFFTScale = 1. / static_cast<double>(iInterpDet * iInterpDet);
760 int iMidPoint = iInterpDet / 2;
761 double dMidPoint = static_cast<double>(iInterpDet) / 2.;
762 int iZerosAdded = iNumInterpDetWithZeros - iInterpDet;
764 // For each view, interpolate to nx length, shift to center at origin, and FFt transform
765 for (unsigned int iView = 0; iView < m_nView; iView++) {
766 DetectorValue* detval = pProj->getDetectorArray(iView).detValues();
767 LinearInterpolator<DetectorValue> projInterp (detval, pProj->m_nDet);
768 for (unsigned int iDet = 0; iDet < iInterpDet; iDet++) {
769 double dInterpPos = (m_nDet / 2.) + (iDet - dMidPoint) * dInterpScale;
770 pcIn[iDet].re = projInterp.interpolate (dInterpPos) * dInterpScale;
774 Fourier::shuffleFourierToNaturalOrder (pcIn, iInterpDet);
775 if (iZerosAdded > 0) {
776 for (unsigned int iDet1 = iMidPoint; iDet1 < iInterpDet; iDet1++)
777 pcIn[iDet1+iZerosAdded] = pcIn[iDet1];
778 for (unsigned int iDet2 = iMidPoint; iDet2 < iMidPoint + iZerosAdded; iDet2++)
779 pcIn[iDet2].re = pcIn[iDet2].im = 0;
782 fftw_one (plan, pcIn, NULL);
784 ppcDetValue[iView] = new std::complex<double> [iNumInterpDetWithZeros];
785 for (unsigned int iD = 0; iD < iNumInterpDetWithZeros; iD++) {
786 ppcDetValue[iView][iD] = std::complex<double> (pcIn[iD].re * dFFTScale, pcIn[iD].im * dFFTScale);
789 Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], iNumInterpDetWithZeros);
793 fftw_destroy_plan (plan);
795 Array2d<double> adView (nx, ny);
796 Array2d<double> adDet (nx, ny);
797 double** ppdView = adView.getArray();
798 double** ppdDet = adDet.getArray();
799 pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet, iNumInterpDetWithZeros, dZeropadRatio,
800 pProj->m_detInc * dInterpScale);
802 pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet,
803 iNumInterpDetWithZeros, iInterpolationID);
805 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
808 for (int i = 0; i < m_nView; i++)
809 delete [] ppcDetValue[i];
810 delete [] ppcDetValue;
818 Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double** ppdView, double** ppdDet,
819 int iNumDetWithZeros, double dZeropadRatio, double dDetInc)
821 // double dLength = viewDiameter();
822 double dLength = phmLen();
823 double xMin = -dLength / 2;
824 double xMax = xMin + dLength;
825 double yMin = -dLength / 2;
826 double yMax = yMin + dLength;
827 double xCent = (xMin + xMax) / 2;
828 double yCent = (yMin + yMax) / 2;
830 xMin = (xMin - xCent) * dZeropadRatio + xCent;
831 xMax = (xMax - xCent) * dZeropadRatio + xCent;
832 yMin = (yMin - yCent) * dZeropadRatio + yCent;
833 yMax = (yMax - yCent) * dZeropadRatio + yCent;
835 double xInc = (xMax - xMin) / nx; // size of cells
836 double yInc = (yMax - yMin) / ny;
838 // +1 is correct for frequency data, ndet-1 is correct for projections
839 int iDetCenter = (iNumDetWithZeros - 1) / 2; // index refering to L=0 projection
840 if (isEven (iNumDetWithZeros))
841 iDetCenter = (iNumDetWithZeros + 1) / 2;
843 // Calculates polar coordinates (view#, det#) for each point on phantom grid
844 double x = xMin + xInc / 2; // Rectang coords of center of pixel
845 for (unsigned int ix = 0; ix < nx; x += xInc, ix++) {
846 double y = yMin + yInc / 2;
847 for (unsigned int iy = 0; iy < ny; y += yInc, iy++) {
848 double r = ::sqrt (x * x + y * y);
849 double phi = atan2 (y, x);
858 ppdView[ix][iy] = (phi - m_rotStart) / m_rotInc;
859 ppdDet[ix][iy] = (r / dDetInc) + iDetCenter;
865 Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag,
866 unsigned int nx, unsigned int ny, std::complex<double>** ppcDetValue, double** ppdView,
867 double** ppdDet, unsigned int nView, unsigned int nDet, unsigned int nDetWithZeros, int iInterpolationID)
869 typedef std::complex<double> complexValue;
871 BilinearInterpolator<complexValue>* pBilinear;
872 if (iInterpolationID == POLAR_INTERP_BILINEAR)
873 pBilinear = new BilinearInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
875 BicubicPolyInterpolator<complexValue>* pBicubic;
876 if (iInterpolationID == POLAR_INTERP_BICUBIC)
877 pBicubic = new BicubicPolyInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
879 for (unsigned int ix = 0; ix < ny; ix++) {
880 for (unsigned int iy = 0; iy < ny; iy++) {
882 if (iInterpolationID == POLAR_INTERP_NEAREST) {
883 unsigned int iView = nearest<int> (ppdView[ix][iy]);
884 unsigned int iDet = nearest<int> (ppdDet[ix][iy]);
885 if (iView == nView) {
887 iDet = m_nDet - iDet;
889 if (iDet >= 0 && iDet < nDetWithZeros && iView >= 0 && iView < nView) {
890 v[ix][iy] = ppcDetValue[iView][iDet].real();
892 vImag[ix][iy] = ppcDetValue[iView][iDet].imag();
896 } else if (iInterpolationID == POLAR_INTERP_BILINEAR) {
897 std::complex<double> vInterp = pBilinear->interpolate (ppdView[ix][iy], ppdDet[ix][iy]);
898 v[ix][iy] = vInterp.real();
900 vImag[ix][iy] = vInterp.imag();
901 } else if (iInterpolationID == POLAR_INTERP_BICUBIC) {
902 std::complex<double> vInterp = pBicubic->interpolate (ppdView[ix][iy], ppdDet[ix][iy]);
903 v[ix][iy] = vInterp.real();
905 vImag[ix][iy] = vInterp.imag();
912 Projections::initFromSomatomAR_STAR (int iNViews, int iNDets, unsigned char* pData, unsigned long lDataLength)
914 init (iNViews, iNDets);
915 m_geometry = Scanner::GEOMETRY_EQUIANGULAR;
916 m_dFocalLength = 510;
917 m_dSourceDetectorLength = 890;
918 m_detInc = convertDegreesToRadians (3.06976 / 60);
919 m_dFanBeamAngle = iNDets * m_detInc;
920 m_detStart = -(m_dFanBeamAngle / 2);
921 m_rotInc = TWOPI / static_cast<double>(iNViews);
923 m_dViewDiameter = sin (m_dFanBeamAngle / 2) * m_dFocalLength * 2;
925 if (! ((iNViews == 750 && lDataLength == 1560000L) || (iNViews == 950 && lDataLength == 1976000L)
926 || (iNViews == 1500 && lDataLength == 3120000)))
929 double dCenter = (iNDets - 1.) / 2.; // change from (Nm+1)/2 because of 0 vs. 1 indexing
930 double* pdCosScale = new double [iNDets];
931 for (int i = 0; i < iNDets; i++)
932 pdCosScale[i] = 1. / cos ((i - dCenter) * m_detInc);
935 for (int iv = 0; iv < iNViews; iv++) {
936 unsigned char* pArgBase = pData + lDataPos;
937 unsigned char* p = pArgBase+0; SwapBytes4IfLittleEndian (p);
938 long lProjNumber = *reinterpret_cast<long*>(p);
940 p = pArgBase+20; SwapBytes4IfLittleEndian (p);
941 long lEscale = *reinterpret_cast<long*>(p);
943 p = pArgBase+28; SwapBytes4IfLittleEndian (p);
944 long lTime = *reinterpret_cast<long*>(p);
946 p = pArgBase + 4; SwapBytes4IfLittleEndian (p);
947 double dAlpha = *reinterpret_cast<float*>(p) + HALFPI;
949 p = pArgBase+12; SwapBytes4IfLittleEndian (p);
950 double dAlign = *reinterpret_cast<float*>(p);
952 p = pArgBase + 16; SwapBytes4IfLittleEndian (p);
953 double dMaxValue = *reinterpret_cast<float*>(p);
955 DetectorArray& detArray = getDetectorArray (iv);
956 detArray.setViewAngle (dAlpha);
957 DetectorValue* detval = detArray.detValues();
959 double dViewScale = 1. / (2294.4871 * ::pow (2.0, -lEscale));
961 for (int id = 0; id < iNDets; id++) {
962 int iV = pData[lDataPos+1] + (pData[lDataPos] << 8);
963 if (iV > 32767) // two's complement signed conversion
965 detval[id] = iV * dViewScale * pdCosScale[id];
969 for (int k = iNDets - 2; k >= 0; k--)
970 detval[k+1] = detval[k];
980 Projections::interpolateToParallel () const
982 if (m_geometry == Scanner::GEOMETRY_PARALLEL)
983 return const_cast<Projections*>(this);
987 Projections* pProjNew = new Projections (nView, nDet);
988 pProjNew->m_geometry = Scanner::GEOMETRY_PARALLEL;
989 pProjNew->m_dFocalLength = m_dFocalLength;
990 pProjNew->m_dSourceDetectorLength = m_dSourceDetectorLength;
991 pProjNew->m_dViewDiameter = m_dViewDiameter;
992 pProjNew->m_dFanBeamAngle = m_dFanBeamAngle;
993 pProjNew->m_calcTime = 0;
994 pProjNew->m_remark = m_remark;
995 pProjNew->m_remark += "; Interpolate to Parallel";
996 pProjNew->m_label.setLabelType (Array2dFileLabel::L_HISTORY);
997 pProjNew->m_label.setLabelString (pProjNew->m_remark);
998 pProjNew->m_label.setCalcTime (pProjNew->m_calcTime);
999 pProjNew->m_label.setDateTime (pProjNew->m_year, pProjNew->m_month, pProjNew->m_day, pProjNew->m_hour, pProjNew->m_minute, pProjNew->m_second);
1001 pProjNew->m_rotStart = 0;
1002 #ifdef CONVERT_PARALLEL_PI
1003 pProjNew->m_rotInc = PI / nView;;
1005 pProjNew->m_rotInc = TWOPI / nView;
1007 pProjNew->m_detStart = -m_dViewDiameter / 2;
1008 pProjNew->m_detInc = m_dViewDiameter / nDet;
1009 if (isEven (nDet)) // even
1010 pProjNew->m_detInc = m_dViewDiameter / (nDet - 1);
1012 ParallelRaysums parallel (this, ParallelRaysums::THETA_RANGE_NORMALIZE_TO_TWOPI);
1014 double* pdThetaValuesForT = new double [pProjNew->nView()];
1015 double* pdRaysumsForT = new double [pProjNew->nView()];
1017 // interpolate to evenly spaced theta (views)
1018 double dDetPos = pProjNew->m_detStart;
1019 for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc) {
1020 parallel.getThetaAndRaysumsForT (iD, pdThetaValuesForT, pdRaysumsForT);
1022 double dViewAngle = m_rotStart;
1023 int iLastFloor = -1;
1024 for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
1025 DetectorValue* detValues = pProjNew->getDetectorArray (iV).detValues();
1026 LinearInterpolator<double> interp (pdThetaValuesForT, pdRaysumsForT, pProjNew->nView(), false);
1027 detValues[iD] = interp.interpolate (dViewAngle, &iLastFloor);
1030 delete pdThetaValuesForT;
1031 delete pdRaysumsForT;
1033 // interpolate to evenly space t (detectors)
1034 double* pdOriginalDetPositions = new double [pProjNew->nDet()];
1035 parallel.getDetPositions (pdOriginalDetPositions);
1037 double* pdDetValueCopy = new double [pProjNew->nDet()];
1038 double dViewAngle = m_rotStart;
1039 for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
1040 DetectorArray& detArray = pProjNew->getDetectorArray (iV);
1041 DetectorValue* detValues = detArray.detValues();
1042 detArray.setViewAngle (dViewAngle);
1044 for (int i = 0; i < pProjNew->nDet(); i++)
1045 pdDetValueCopy[i] = detValues[i];
1047 double dDetPos = pProjNew->m_detStart;
1048 int iLastFloor = -1;
1049 LinearInterpolator<double> interp (pdOriginalDetPositions, pdDetValueCopy, pProjNew->nDet(), false);
1050 for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc)
1051 detValues[iD] = interp.interpolate (dDetPos, &iLastFloor);
1053 delete pdDetValueCopy;
1054 delete pdOriginalDetPositions;
1060 ///////////////////////////////////////////////////////////////////////////////
1062 // Class ParallelRaysums
1064 // Used for converting divergent beam raysums into Parallel raysums
1066 ///////////////////////////////////////////////////////////////////////////////
1068 ParallelRaysums::ParallelRaysums (const Projections* pProjections, int iThetaRange)
1069 : m_iNumCoordinates(0), m_iNumView(pProjections->nView()), m_iNumDet(pProjections->nDet()),
1070 m_iThetaRange (iThetaRange), m_pCoordinates(NULL)
1072 int iGeometry = pProjections->geometry();
1073 double dDetInc = pProjections->detInc();
1074 double dDetStart = pProjections->detStart();
1075 double dFocalLength = pProjections->focalLength();
1077 m_iNumCoordinates = m_iNumView * m_iNumDet;
1078 m_pCoordinates = new ParallelRaysumCoordinate [m_iNumCoordinates];
1079 m_vecpCoordinates.reserve (m_iNumCoordinates);
1080 for (int i = 0; i < m_iNumCoordinates; i++)
1081 m_vecpCoordinates[i] = m_pCoordinates + i;
1083 int iCoordinate = 0;
1084 for (int iV = 0; iV < m_iNumView; iV++) {
1085 double dViewAngle = pProjections->getDetectorArray(iV).viewAngle();
1086 const DetectorValue* detValues = pProjections->getDetectorArray(iV).detValues();
1088 double dDetPos = dDetStart;
1089 for (int iD = 0; iD < m_iNumDet; iD++) {
1090 ParallelRaysumCoordinate* pC = m_vecpCoordinates[iCoordinate++];
1092 if (iGeometry == Scanner::GEOMETRY_PARALLEL) {
1093 pC->m_dTheta = dViewAngle;
1095 } else if (iGeometry == Scanner::GEOMETRY_EQUILINEAR) {
1096 double dFanAngle = atan (dDetPos / pProjections->sourceDetectorLength());
1097 pC->m_dTheta = dViewAngle + dFanAngle;
1098 pC->m_dT = dFocalLength * sin(dFanAngle);
1100 } else if (iGeometry == Scanner::GEOMETRY_EQUIANGULAR) {
1101 // fan angle is same as dDetPos
1102 pC->m_dTheta = dViewAngle + dDetPos;
1103 pC->m_dT = dFocalLength * sin (dDetPos);
1105 if (m_iThetaRange != THETA_RANGE_UNCONSTRAINED) {
1106 pC->m_dTheta = normalizeAngle (pC->m_dTheta);
1107 if (m_iThetaRange == THETA_RANGE_FOLD_TO_PI && pC->m_dTheta >= PI) {
1109 pC->m_dT = -pC->m_dT;
1112 pC->m_dRaysum = detValues[iD];
1118 ParallelRaysums::~ParallelRaysums()
1120 delete m_pCoordinates;
1123 ParallelRaysums::CoordinateContainer&
1124 ParallelRaysums::getSortedByTheta()
1126 if (m_vecpSortedByTheta.size() == 0) {
1127 m_vecpSortedByTheta.resize (m_iNumCoordinates);
1128 for (int i = 0; i < m_iNumCoordinates; i++)
1129 m_vecpSortedByTheta[i] = m_vecpCoordinates[i];
1130 std::sort (m_vecpSortedByTheta.begin(), m_vecpSortedByTheta.end(), ParallelRaysumCoordinate::compareByTheta);
1133 return m_vecpSortedByTheta;
1136 ParallelRaysums::CoordinateContainer&
1137 ParallelRaysums::getSortedByT()
1139 if (m_vecpSortedByT.size() == 0) {
1140 m_vecpSortedByT.resize (m_iNumCoordinates);
1141 for (int i = 0; i < m_iNumCoordinates; i++)
1142 m_vecpSortedByT[i] = m_vecpCoordinates[i];
1143 std::sort (m_vecpSortedByT.begin(), m_vecpSortedByT.end(), ParallelRaysumCoordinate::compareByT);
1146 return m_vecpSortedByT;
1151 ParallelRaysums::getLimits (double* dMinT, double* dMaxT, double* dMinTheta, double* dMaxTheta) const
1153 if (m_iNumCoordinates <= 0)
1156 *dMinT = *dMaxT = m_vecpCoordinates[0]->m_dT;
1157 *dMinTheta = *dMaxTheta = m_vecpCoordinates[0]->m_dTheta;
1159 for (int i = 0; i < m_iNumCoordinates; i++) {
1160 double dT = m_vecpCoordinates[i]->m_dT;
1161 double dTheta = m_vecpCoordinates[i]->m_dTheta;
1165 else if (dT > *dMaxT)
1168 if (dTheta < *dMinTheta)
1169 *dMinTheta = dTheta;
1170 else if (dTheta > *dMaxTheta)
1171 *dMaxTheta = dTheta;
1176 ParallelRaysums::getThetaAndRaysumsForT (int iTheta, double* pTheta, double* pRaysum)
1178 const CoordinateContainer& coordsT = getSortedByT();
1180 int iBase = iTheta * m_iNumView;
1181 for (int i = 0; i < m_iNumView; i++) {
1182 int iPos = iBase + i;
1183 pTheta[i] = coordsT[iPos]->m_dTheta;
1184 pRaysum[i] = coordsT[iPos]->m_dRaysum;
1189 ParallelRaysums::getDetPositions (double* pdDetPos)
1191 const CoordinateContainer& coordsT = getSortedByT();
1194 for (int i = 0; i < m_iNumDet; i++) {
1195 pdDetPos[i] = coordsT[iPos]->m_dT;