r636: Optimized Rebinning, Added Reconstruct with Rebinning option
[ctsim.git] / libctsim / projections.cpp
1 /*****************************************************************************
2 ** FILE IDENTIFICATION
3 **
4 **   Name:         projections.cpp         Projection data classes
5 **   Programmer:   Kevin Rosenberg
6 **   Date Started: Aug 84
7 **
8 **  This is part of the CTSim program
9 **  Copyright (c) 1983-2001 Kevin Rosenberg
10 **
11 **  $Id: projections.cpp,v 1.63 2001/03/13 08:24:41 kevin Exp $
12 **
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.
16 **
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.
21 **
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 ******************************************************************************/
26
27 #include "ct.h"
28
29 const kuint16 Projections::m_signature = ('P'*256 + 'J');
30
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;
35
36 const char* const Projections::s_aszInterpName[] = 
37 {
38   {"nearest"},
39   {"bilinear"},
40 //  {"bicubic"},
41 };
42
43 const char* const Projections::s_aszInterpTitle[] = 
44 {
45   {"Nearest"},
46   {"Bilinear"},
47 //  {"Bicubic"},
48 };
49
50 const int Projections::s_iInterpCount = sizeof(s_aszInterpName) / sizeof(char*);
51
52
53 /* NAME
54 *    Projections                Constructor for projections matrix storage 
55 *
56 * SYNOPSIS
57 *    proj = projections_create (filename, nView, nDet)
58 *    Projections& proj          Allocated projections structure & matrix
59 *    int nView                  Number of rotated view
60 *    int nDet                   Number of detectors
61 *
62 */
63
64 Projections::Projections (const Scanner& scanner)
65 : m_projData(0)
66 {
67   initFromScanner (scanner);
68 }
69
70
71 Projections::Projections (const int nView, const int nDet)
72 : m_projData(0)
73 {
74   init (nView, nDet);
75 }
76
77 Projections::Projections (void)
78 : m_projData(0)
79 {
80   init (0, 0);
81 }
82
83 Projections::~Projections (void)
84 {
85   deleteProjData();
86 }
87
88 int
89 Projections::convertInterpNameToID (const char* const interpName)
90 {
91   int interpID = POLAR_INTERP_INVALID;
92   
93   for (int i = 0; i < s_iInterpCount; i++)
94     if (strcasecmp (interpName, s_aszInterpName[i]) == 0) {
95       interpID = i;
96       break;
97     }
98     
99     return (interpID);
100 }
101
102 const char*
103 Projections::convertInterpIDToName (const int interpID)
104 {
105   static const char *interpName = "";
106   
107   if (interpID >= 0 && interpID < s_iInterpCount)
108     return (s_aszInterpName[interpID]);
109   
110   return (interpName);
111 }
112
113 const char*
114 Projections::convertInterpIDToTitle (const int interpID)
115 {
116   static const char *interpTitle = "";
117   
118   if (interpID >= 0 && interpID < s_iInterpCount)
119     return (s_aszInterpTitle[interpID]);
120   
121   return (interpTitle);
122 }
123
124
125
126 void
127 Projections::init (const int nView, const int nDet)
128 {
129   m_label.setLabelType (Array2dFileLabel::L_HISTORY);
130   m_nView = nView;
131   m_nDet = nDet;
132   newProjData ();
133   
134   time_t t = time (NULL);
135   tm* lt = localtime (&t);
136   m_year = lt->tm_year;
137   m_month = lt->tm_mon;
138   m_day = lt->tm_mday;
139   m_hour = lt->tm_hour;
140   m_minute = lt->tm_min;
141   m_second = lt->tm_sec;
142 }
143
144 void
145 Projections::initFromScanner (const Scanner& scanner)
146 {
147   m_label.setLabelType (Array2dFileLabel::L_HISTORY);
148   deleteProjData();
149   init (scanner.nView(), scanner.nDet());
150   
151   m_rotInc = scanner.rotInc();
152   m_detInc = scanner.detInc();
153   m_detStart =  scanner.detStart();
154   m_geometry = scanner.geometry();
155   m_dFocalLength = scanner.focalLength();
156   m_dSourceDetectorLength = scanner.sourceDetectorLength();
157   m_dViewDiameter = scanner.viewDiameter();
158   m_rotStart = 0;
159   m_dFanBeamAngle = scanner.fanBeamAngle();
160 }
161
162 void
163 Projections::setNView (int nView)  // used by MPI to reduce # of views
164 {
165   deleteProjData();
166   init (nView, m_nDet);
167 }
168
169 // NAME
170 // newProjData
171
172 void 
173 Projections::newProjData (void)
174 {
175   if (m_projData)
176     sys_error(ERR_WARNING, "m_projData != NULL [newProjData]");
177   
178   if (m_nView > 0 && m_nDet) {
179     m_projData = new DetectorArray* [m_nView];
180     
181     for (int i = 0; i < m_nView; i++)
182       m_projData[i] = new DetectorArray (m_nDet);
183   }
184 }
185
186
187 /* NAME
188 *   projections_free                    Free memory allocated to projections
189 *
190 * SYNOPSIS
191 *   projections_free(proj)
192 *   Projections& proj                           Projectionss to be deallocated
193 */
194
195 void 
196 Projections::deleteProjData (void)
197 {
198   if (m_projData != NULL) {
199     for (int i = 0; i < m_nView; i++)
200       delete m_projData[i];
201     
202     delete m_projData;
203     m_projData = NULL;
204   }
205 }
206
207
208 /* NAME
209 *    Projections::headerWwrite         Write data header for projections file
210 *
211 */
212
213 bool 
214 Projections::headerWrite (fnetorderstream& fs)
215 {
216   kuint16 _hsize = m_headerSize;
217   kuint16 _signature = m_signature;
218   kuint32 _nView = m_nView;
219   kuint32 _nDet = m_nDet;
220   kuint32 _geom = m_geometry;
221   kuint16 _remarksize = m_remark.length();
222   kuint16 _year = m_year;
223   kuint16 _month = m_month;
224   kuint16 _day = m_day;
225   kuint16 _hour = m_hour;
226   kuint16 _minute = m_minute;
227   kuint16 _second = m_second;
228   
229   kfloat64 _calcTime = m_calcTime;
230   kfloat64 _rotStart = m_rotStart;
231   kfloat64 _rotInc = m_rotInc;
232   kfloat64 _detStart = m_detStart;
233   kfloat64 _detInc = m_detInc;
234   kfloat64 _viewDiameter = m_dViewDiameter;
235   kfloat64 _focalLength = m_dFocalLength;
236   kfloat64 _sourceDetectorLength = m_dSourceDetectorLength;
237   kfloat64 _fanBeamAngle = m_dFanBeamAngle;
238
239   fs.seekp(0);
240   if (! fs)
241     return false;
242   
243   fs.writeInt16 (_hsize);
244   fs.writeInt16 (_signature);
245   fs.writeInt32 (_nView);
246   fs.writeInt32 (_nDet);
247   fs.writeInt32 (_geom);
248   fs.writeFloat64 (_calcTime);
249   fs.writeFloat64 (_rotStart);
250   fs.writeFloat64 (_rotInc);
251   fs.writeFloat64 (_detStart);
252   fs.writeFloat64 (_detInc);
253   fs.writeFloat64 (_viewDiameter);
254   fs.writeFloat64 (_focalLength);
255   fs.writeFloat64 (_sourceDetectorLength);
256   fs.writeFloat64 (_fanBeamAngle);
257   fs.writeInt16 (_year);
258   fs.writeInt16 (_month);
259   fs.writeInt16 (_day);
260   fs.writeInt16 (_hour);
261   fs.writeInt16 (_minute);
262   fs.writeInt16 (_second);
263   fs.writeInt16 (_remarksize);
264   fs.write (m_remark.c_str(), _remarksize);
265   
266   m_headerSize = fs.tellp();
267   _hsize = m_headerSize;
268   fs.seekp(0);
269   fs.writeInt16 (_hsize);
270   if (! fs)
271     return false;
272   
273   return true;
274 }
275
276 /* NAME
277 *    projections_read_header         Read data header for projections file
278 *
279 */
280 bool
281 Projections::headerRead (fnetorderstream& fs)
282 {
283   kuint16 _hsize, _signature, _year, _month, _day, _hour, _minute, _second, _remarksize = 0;
284   kuint32 _nView, _nDet, _geom;
285   kfloat64 _calcTime, _rotStart, _rotInc, _detStart, _detInc, _focalLength, _sourceDetectorLength, _viewDiameter, _fanBeamAngle;
286   
287   fs.seekg(0);
288   if (! fs)
289     return false;
290   
291   fs.readInt16 (_hsize);
292   fs.readInt16 (_signature);
293   fs.readInt32 (_nView);
294   fs.readInt32 (_nDet);
295   fs.readInt32 (_geom);
296   fs.readFloat64 (_calcTime);
297   fs.readFloat64 (_rotStart);
298   fs.readFloat64 (_rotInc);
299   fs.readFloat64 (_detStart);
300   fs.readFloat64 (_detInc);
301   fs.readFloat64 (_viewDiameter);
302   fs.readFloat64 (_focalLength);
303   fs.readFloat64 (_sourceDetectorLength);
304   fs.readFloat64 (_fanBeamAngle);
305   fs.readInt16 (_year);
306   fs.readInt16 (_month);
307   fs.readInt16 (_day);
308   fs.readInt16 (_hour);
309   fs.readInt16 (_minute);
310   fs.readInt16 (_second);
311   fs.readInt16 (_remarksize);
312   
313   if (! fs) {
314     sys_error (ERR_SEVERE, "Error reading header information , _remarksize=%d [projections_read_header]", _remarksize);
315     return false;
316   }
317   
318   if (_signature != m_signature) {
319     sys_error (ERR_SEVERE, "File %s does not have a valid projection file signature", m_filename.c_str());
320     return false;
321   }
322   
323   char* pszRemarkStorage = new char [_remarksize+1];
324   fs.read (pszRemarkStorage, _remarksize);
325   if (! fs) {
326     sys_error (ERR_SEVERE, "Error reading remark, _remarksize = %d", _remarksize);
327     return false;
328   }
329   pszRemarkStorage[_remarksize] = 0;
330   m_remark = pszRemarkStorage;
331   delete pszRemarkStorage;
332   
333   off_t _hsizeread = fs.tellg();
334   if (!fs || _hsizeread != _hsize) {
335     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);
336     return false;
337   }
338   
339   m_headerSize = _hsize;
340   m_nView = _nView;
341   m_nDet = _nDet;
342   m_geometry = _geom;
343   m_calcTime = _calcTime;
344   m_rotStart = _rotStart;
345   m_rotInc = _rotInc;
346   m_detStart = _detStart;
347   m_detInc = _detInc;
348   m_dFocalLength = _focalLength;
349   m_dSourceDetectorLength = _sourceDetectorLength;
350   m_dViewDiameter = _viewDiameter;
351   m_dFanBeamAngle = _fanBeamAngle;
352   m_year = _year;
353   m_month = _month;
354   m_day = _day;
355   m_hour = _hour;
356   m_minute = _minute;
357   m_second = _second;
358   
359   m_label.setLabelType (Array2dFileLabel::L_HISTORY);
360   m_label.setLabelString (m_remark);
361   m_label.setCalcTime (m_calcTime);
362   m_label.setDateTime (m_year, m_month, m_day, m_hour, m_minute, m_second);
363   
364   return true;
365 }
366
367 bool
368 Projections::read (const std::string& filename)
369 {
370   return read (filename.c_str());
371 }
372
373
374 bool
375 Projections::read (const char* filename)
376 {
377   m_filename = filename;
378 #ifdef MSVC
379   frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary);
380 #else
381   frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary | std::ios::nocreate);
382 #endif
383   
384   if (fileRead.fail())
385     return false;
386   
387   if (! headerRead (fileRead))
388     return false;
389   
390   deleteProjData ();
391   newProjData();
392   
393   for (int i = 0; i < m_nView; i++) {
394     if (! detarrayRead (fileRead, *m_projData[i], i))
395       break;
396   }
397   
398   fileRead.close();
399   return true;
400 }
401
402
403 bool 
404 Projections::copyViewData (const std::string& filename, std::ostream& os, int startView, int endView)
405 {
406   return copyViewData (filename.c_str(), os, startView, endView);
407 }
408
409 bool 
410 Projections::copyViewData (const char* const filename, std::ostream& os, int startView, int endView)
411 {
412   frnetorderstream is (filename, std::ios::in | std::ios::binary);
413   kuint16 sizeHeader, signature;
414   kuint32 _nView, _nDet;
415   
416   is.seekg (0);
417   if (is.fail()) {
418     sys_error (ERR_SEVERE, "Unable to read projection file %s", filename);
419     return false;
420   }
421
422   is.readInt16 (sizeHeader);
423   is.readInt16 (signature);
424   is.readInt32 (_nView);
425   is.readInt32 (_nDet);
426   int nView = _nView;
427   int nDet = _nDet;
428   
429   if (signature != m_signature) {
430     sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
431     return false;
432   }
433   
434   if (startView < 0)
435     startView = 0;
436   if (startView > nView - 1)
437     startView = nView;
438   if (endView < 0 || endView > nView - 1)
439     endView = nView - 1;
440   
441   if (startView > endView) { // swap if start > end
442     int tempView = endView;
443     endView = startView;
444     startView = tempView;
445   }
446   
447   int sizeView = 8 /* view_angle */ + 4 /* nDet */ + (4 * nDet);
448   unsigned char* pViewData = new unsigned char [sizeView];
449   
450   for (int i = startView; i <= endView; i++) {
451     is.seekg (sizeHeader + i * sizeView);
452     is.read (reinterpret_cast<char*>(pViewData), sizeView);
453     os.write (reinterpret_cast<char*>(pViewData), sizeView);
454     if (is.fail() || os.fail())
455       break;
456   }
457   
458   delete pViewData;
459   if (is.fail()) 
460     sys_error (ERR_SEVERE, "Error reading projection file");
461   if (os.fail()) 
462     sys_error (ERR_SEVERE, "Error writing projection file");
463   
464   return (! (is.fail() | os.fail()));
465 }
466
467 bool 
468 Projections::copyHeader (const std::string& filename, std::ostream& os)
469 {
470   return copyHeader (filename.c_str(), os);
471 }
472
473 bool
474 Projections::copyHeader (const char* const filename, std::ostream& os)
475 {
476   frnetorderstream is (filename, std::ios::in | std::ios::binary);
477   kuint16 sizeHeader, signature;
478   is.readInt16 (sizeHeader);
479   is.readInt16 (signature);
480   is.seekg (0);
481   if (signature != m_signature) {
482     sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
483     return false;
484   }
485   
486   unsigned char* pHdrData = new unsigned char [sizeHeader];
487   is.read (reinterpret_cast<char*>(pHdrData), sizeHeader);
488   if (is.fail()) {
489     sys_error (ERR_SEVERE, "Error reading header");
490     return false;
491   }
492   
493   os.write (reinterpret_cast<char*>(pHdrData), sizeHeader);
494   if (os.fail()) {
495     sys_error (ERR_SEVERE, "Error writing header");
496     return false;
497   }
498   
499   return true;
500 }
501
502 bool
503 Projections::write (const std::string& filename)
504 {
505   return write (filename.c_str());
506 }
507
508 bool
509 Projections::write (const char* filename)
510 {
511   frnetorderstream fs (filename, std::ios::out | std::ios::binary | std::ios::trunc | std::ios::ate);
512   m_filename = filename;
513   if (! fs) {
514     sys_error (ERR_SEVERE, "Error opening file %s for output [projections_create]", filename);
515     return false;
516   }
517   
518   if (! headerWrite (fs))
519     return false;
520   
521   if (m_projData != NULL) {
522     for (int i = 0; i < m_nView; i++) {
523       if (! detarrayWrite (fs, *m_projData[i], i))
524         break;
525     }
526   }
527   if (! fs)
528     return false;
529   
530   fs.close();
531   
532   return true;
533 }
534
535 /* NAME
536 *   detarrayRead                Read a Detector Array structure from the disk
537 *
538 * SYNOPSIS
539 *   detarrayRead (proj, darray, view_num)
540 *   DETARRAY *darray            Detector array storage location to be filled
541 *   int      view_num           View number to read
542 */
543
544 bool
545 Projections::detarrayRead (fnetorderstream& fs, DetectorArray& darray, const int iview)
546 {
547   const int detval_bytes = darray.nDet() * sizeof(kfloat32);
548   const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
549   const int view_bytes = detheader_bytes + detval_bytes;
550   const off_t start_data = m_headerSize + (iview * view_bytes);
551   DetectorValue* detval_ptr = darray.detValues();  
552   kfloat64 view_angle;
553   kuint32 nDet;
554   
555   fs.seekg (start_data);
556   
557   fs.readFloat64 (view_angle);
558   fs.readInt32 (nDet);
559   darray.setViewAngle (view_angle);
560   //  darray.setNDet ( nDet);
561   
562   for (unsigned int i = 0; i < nDet; i++) {
563     kfloat32 detval;
564     fs.readFloat32 (detval);
565     detval_ptr[i] = detval;
566   }
567   if (! fs)
568     return false;
569   
570   return true;
571 }
572
573
574 /* NAME
575 *   detarrayWrite                       Write detector array data to the disk
576 *
577 * SYNOPSIS
578 *   detarrayWrite (darray, view_num)
579 *   DETARRAY *darray                    Detector array data to be written
580 *   int      view_num                   View number to write
581 *
582 * DESCRIPTION
583 *       This routine writes the detarray data from the disk sequentially to
584 *    the file that was opened with open_projections().  Data is written in
585 *    binary format.
586 */
587
588 bool
589 Projections::detarrayWrite (fnetorderstream& fs, const DetectorArray& darray, const int iview)
590 {
591   const int detval_bytes = darray.nDet() * sizeof(float);
592   const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
593   const int view_bytes = detheader_bytes + detval_bytes;
594   const off_t start_data = m_headerSize + (iview * view_bytes);
595   const DetectorValue* const detval_ptr = darray.detValues();  
596   kfloat64 view_angle = darray.viewAngle();
597   kuint32 nDet = darray.nDet();
598   
599   fs.seekp (start_data);
600   if (! fs) {
601     sys_error (ERR_SEVERE, "Error seeking detectory array [detarrayWrite]");
602     return false;
603   }
604   
605   fs.writeFloat64 (view_angle);
606   fs.writeInt32 (nDet);
607   
608   for (unsigned int i = 0; i < nDet; i++) {
609     kfloat32 detval = detval_ptr[i];
610     fs.writeFloat32 (detval);
611   }
612   
613   if (! fs)
614     return (false);
615   
616   return true;
617 }
618
619 /* NAME
620 *   printProjectionData                 Print projections data
621 *
622 * SYNOPSIS
623 *   printProjectionData ()
624 */
625
626 void
627 Projections::printProjectionData ()
628 {
629   printProjectionData (0, nView() - 1);
630 }
631
632 void
633 Projections::printProjectionData (int startView, int endView)
634 {
635   printf("Projections Data\n\n");
636   printf("Description: %s\n", m_remark.c_str());
637   printf("Geometry: %s\n", Scanner::convertGeometryIDToName (m_geometry));
638   printf("nView       = %8d             nDet = %8d\n", m_nView, m_nDet);
639   printf("focalLength = %8.4f   ViewDiameter = %8.4f\n", m_dFocalLength, m_dViewDiameter);
640   printf("fanBeamAngle= %8.4f SourceDetector = %8.4f\n", convertRadiansToDegrees(m_dFanBeamAngle), m_dSourceDetectorLength);
641   printf("rotStart    = %8.4f         rotInc = %8.4f\n", m_rotStart, m_rotInc);
642   printf("detStart    = %8.4f         detInc = %8.4f\n", m_detStart, m_detInc);
643   if (m_projData != NULL) {
644     if (startView < 0)
645       startView = 0;
646     if (endView < 0)
647       endView = m_nView - 1;
648     if (startView > m_nView - 1)
649       startView = m_nView - 1;
650     if (endView > m_nView - 1)
651       endView = m_nView - 1;
652     for (int ir = startView; ir <= endView - 1; ir++) {
653       printf("View %d: angle %f\n", ir, m_projData[ir]->viewAngle());
654       DetectorValue* detval = m_projData[ir]->detValues();
655       for (int id = 0; id < m_projData[ir]->nDet(); id++)
656         printf("%8.4f  ", detval[id]);
657       printf("\n");
658     }
659   }
660 }
661
662 void 
663 Projections::printScanInfo (std::ostringstream& os) const
664 {
665   os << "Number of detectors: " << m_nDet << "\n";
666   os << "Number of views: " << m_nView<< "\n";
667   os << "Description: " << m_remark.c_str()<< "\n";
668   os << "Geometry: " << Scanner::convertGeometryIDToName (m_geometry)<< "\n";
669   os << "Focal Length: " << m_dFocalLength<< "\n";
670   os << "Source Detector Length: " << m_dSourceDetectorLength << "\n";
671   os << "View Diameter: " << m_dViewDiameter<< "\n";
672   os << "Fan Beam Angle: " << convertRadiansToDegrees(m_dFanBeamAngle) << "\n";
673   os << "detStart: " << m_detStart<< "\n";
674   os << "detInc: " << m_detInc<< "\n";
675   os << "rotStart: " << m_rotStart<< "\n";
676   os << "rotInc: " << m_rotInc<< "\n";
677 }
678
679
680 bool 
681 Projections::convertPolar (ImageFile& rIF, int iInterpolationID)
682 {
683   unsigned int nx = rIF.nx();
684   unsigned int ny = rIF.ny();
685   ImageFileArray v = rIF.getArray();
686   ImageFileArray vImag = rIF.getImaginaryArray();
687
688   if (! v || nx == 0 || ny == 0)
689     return false;
690
691   Projections* pProj = this;
692   if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
693     pProj = interpolateToParallel();
694   
695   Array2d<double> adView (nx, ny);
696   Array2d<double> adDet (nx, ny);
697   double** ppdView = adView.getArray();
698   double** ppdDet = adDet.getArray();
699
700   if (! pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet)) 
701     return false;
702
703   std::complex<double>** ppcDetValue = new std::complex<double>* [m_nView];
704   unsigned int iView;
705   for (iView = 0; iView < m_nView; iView++) {
706     ppcDetValue[iView] = new std::complex<double> [m_nDet];
707     for (unsigned int iDet = 0; iDet < m_nDet; iDet++)
708       ppcDetValue[iView][iDet] = std::complex<double>(pProj->getDetectorArray (iView).detValues()[iDet], 0);
709   }
710
711   pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet, iInterpolationID);
712
713   for (iView = 0; iView < m_nView; iView++)
714     delete [] ppcDetValue[iView];
715   delete [] ppcDetValue;
716
717   if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
718     delete pProj;
719
720   return true;
721 }
722
723
724 bool 
725 Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad)
726 {
727   unsigned int nx = rIF.nx();
728   unsigned int ny = rIF.ny();
729   ImageFileArray v = rIF.getArray();
730   if (! rIF.isComplex())
731     rIF.convertRealToComplex();
732   ImageFileArray vImag = rIF.getImaginaryArray();
733
734   if (! v || nx == 0 || ny == 0)
735     return false;
736   
737   if (m_geometry != Scanner::GEOMETRY_PARALLEL) {
738     sys_error (ERR_WARNING, "convertFFTPolar supports Parallel only");
739     return false;
740   }
741   
742 #ifndef HAVE_FFT
743   return false;
744 #else
745   Array2d<double> adView (nx, ny);
746   Array2d<double> adDet (nx, ny);
747   double** ppdView = adView.getArray();
748   double** ppdDet = adDet.getArray();
749
750   std::complex<double>** ppcDetValue = new std::complex<double>* [m_nView];
751   unsigned int iView;
752   double* pdDet = new double [m_nDet];
753   fftw_complex* pcIn = new fftw_complex [m_nDet];
754   fftw_plan plan = fftw_create_plan (m_nDet, FFTW_FORWARD, FFTW_IN_PLACE);
755
756   for (iView = 0; iView < m_nView; iView++) {
757     unsigned int iDet;
758     for (iDet = 0; iDet < m_nDet; iDet++) {
759       pcIn[iDet].re = getDetectorArray(iView).detValues()[iDet];
760       pcIn[iDet].im = 0;
761     }
762     fftw_one (plan, pcIn, NULL);
763     ppcDetValue[iView] = new std::complex<double> [m_nDet];
764     for (iDet = 0; iDet < m_nDet; iDet++)
765       ppcDetValue[iView][iDet] = std::complex<double> (pcIn[iDet].re, pcIn[iDet].im); 
766     Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], m_nDet);
767   }
768
769   fftw_destroy_plan (plan);  
770   delete [] pcIn;
771   
772   bool bError = calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet);
773
774   if (! bError)
775     interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, m_nView, m_nDet, iInterpolationID);
776
777   for (iView = 0; iView < m_nView; iView++)
778     delete [] ppcDetValue[iView];
779   delete [] ppcDetValue;
780
781   return bError;
782 #endif
783 }
784
785
786 bool
787 Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double** ppdView, double** ppdDet)
788 {
789   double xMin = -phmLen() / 2;
790   double xMax = xMin + phmLen();
791   double yMin = -phmLen() / 2;
792   double yMax = yMin + phmLen();
793   
794   double xInc = (xMax - xMin) / nx;     // size of cells
795   double yInc = (yMax - yMin) / ny;
796   
797   int iDetCenter = (m_nDet - 1) / 2;    // index refering to L=0 projection 
798
799   // Calculates polar coordinates (view#, det#) for each point on phantom grid
800   double x = xMin + xInc / 2;   // Rectang coords of center of pixel 
801   for (unsigned int ix = 0; ix < nx; x += xInc, ix++) {
802     double y = yMin + yInc / 2;
803     for (unsigned int iy = 0; iy < ny; y += yInc, iy++) {
804       double r = ::sqrt (x * x + y * y);
805       double phi = atan2 (y, x);
806
807       if (phi >= PI) {
808         phi -= PI;
809       } else if (phi < 0) {
810         phi += PI;
811       } else
812         r = -r;
813       
814       ppdView[ix][iy] = (phi - m_rotStart) / m_rotInc;
815       ppdDet[ix][iy] = (r / m_detInc) + iDetCenter;
816     }
817   }
818
819   return true;
820 }
821
822 void
823 Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag,
824      unsigned int nx, unsigned int ny, std::complex<double>** ppcDetValue,
825      double** ppdView, double** ppdDet, unsigned int nView, unsigned int nDet, int iInterpolationID)
826 {
827   for (unsigned int ix = 0; ix < ny; ix++) {
828     for (unsigned int iy = 0; iy < ny; iy++) {
829       if (iInterpolationID == POLAR_INTERP_NEAREST) {
830         unsigned int iView = nearest<int> (ppdView[ix][iy]);
831         unsigned int iDet = nearest<int> (ppdDet[ix][iy]);
832         if (iView == nView) {
833           iView = 0;
834        //   iDet = m_nDet - iDet;
835         }
836         if (iDet >= 0 && iDet < nDet && iView >= 0 && iView < nView) {
837           v[ix][iy] = ppcDetValue[iView][iDet].real();
838           if (vImag)
839             vImag[ix][iy] = ppcDetValue[iView][iDet].imag();
840         } else {
841           sys_error (ERR_SEVERE, "Can't find projection data for ix=%d,iy=%d with radView=%f and radDet=%f", 
842             ix, iy, ppdView[ix][iy], ppdDet[ix][iy]);
843           v[ix][iy] = 0;
844         }
845       } else if (iInterpolationID == POLAR_INTERP_BILINEAR) {
846         unsigned int iFloorView = static_cast<int>(ppdView[ix][iy]);
847         double dFracView = ppdView[ix][iy] - iFloorView;
848         unsigned int iFloorDet = static_cast<int>(ppdDet[ix][iy]);
849         double dFracDet = ppdDet[ix][iy] - iFloorDet;
850
851         if (iFloorDet >= 0 && iFloorView >= 0) { 
852           std::complex<double> v1 = ppcDetValue[iFloorView][iFloorDet];
853           std::complex<double> v2, v3, v4;
854           if (iFloorView < nView - 1)
855             v2 = ppcDetValue[iFloorView + 1][iFloorDet];
856           else 
857             v2 = ppcDetValue[0][iFloorDet];
858           if (iFloorDet < nDet - 1) 
859             v4 = ppcDetValue[iFloorView][iFloorDet+1];
860           else
861             v4 = v1;
862           if (iFloorView < nView - 1 && iFloorDet < nDet - 1)
863             v3 = ppcDetValue [iFloorView+1][iFloorDet+1];
864           else if (iFloorView < nView - 1)
865             v3 = v2;
866           else 
867             v3 = ppcDetValue[0][iFloorDet+1];
868           std::complex<double> vInterp = (1 - dFracView) * (1 - dFracDet) * v1 +
869             dFracView * (1 - dFracDet) * v2 + dFracView * dFracDet * v3 +
870             dFracDet * (1 - dFracView) * v4;
871           v[ix][iy] = vInterp.real();
872           if (vImag)
873             vImag[ix][iy] = vInterp.imag();
874         } else {
875           sys_error (ERR_SEVERE, "Can't find projection data for ix=%d,iy=%d with radView=%f and radDet=%f", 
876             ix, iy, ppdView[ix][iy], ppdDet[ix][iy]);
877           v[ix][iy] = 0;
878           if (vImag)
879             vImag[ix][iy] = 0;
880         }
881       } else if (iInterpolationID == POLAR_INTERP_BICUBIC) {
882         v[ix][iy] =0;
883           if (vImag)
884             vImag[ix][iy] = 0;
885       }
886     }
887   }
888 }
889
890 bool
891 Projections::initFromSomatomAR_STAR (int iNViews, int iNDets, unsigned char* pData, unsigned long lDataLength)
892 {
893   init (iNViews, iNDets);
894   m_geometry = Scanner::GEOMETRY_EQUIANGULAR;
895   m_dFocalLength = 510;
896   m_dSourceDetectorLength = 890;
897   m_detInc = convertDegreesToRadians (3.06976 / 60);
898   m_dFanBeamAngle = (iNDets + 1) * m_detInc;
899   m_detStart = -(m_dFanBeamAngle / 2);
900   m_rotInc = TWOPI / static_cast<double>(iNViews);
901   m_rotStart = HALFPI;
902   m_dViewDiameter = sin (m_dFanBeamAngle / 2) * m_dFocalLength * 2;
903
904   if (! ((iNViews == 750 && lDataLength == 1560000L) || (iNViews == 950 && lDataLength == 1976000L) 
905                 || (iNViews == 1500 && lDataLength == 3120000)))
906     return false;
907
908   double dCenter = (iNDets - 1.) / 2.; // change from (Nm+1)/2 because of 0 vs. 1 indexing
909   double* pdCosScale = new double [iNDets];
910   for (int i = 0; i < iNDets; i++)
911     pdCosScale[i] = 1. / cos ((i - dCenter) * m_detInc);
912
913   long lDataPos = 0;
914   for (int iv = 0; iv < iNViews; iv++) {
915     unsigned char* pArgBase = pData + lDataPos;
916     unsigned char* p = pArgBase+0; SwapBytes4IfLittleEndian (p);
917     long lProjNumber = *reinterpret_cast<long*>(p);
918
919     p = pArgBase+20;  SwapBytes4IfLittleEndian (p);
920     long lEscale = *reinterpret_cast<long*>(p);
921
922     p = pArgBase+28;  SwapBytes4IfLittleEndian (p);
923     long lTime = *reinterpret_cast<long*>(p);
924
925     p = pArgBase + 4; SwapBytes4IfLittleEndian (p);
926     double dAlpha = *reinterpret_cast<float*>(p) + HALFPI;
927
928     p = pArgBase+12; SwapBytes4IfLittleEndian (p);
929     double dAlign = *reinterpret_cast<float*>(p);
930
931     p = pArgBase + 16; SwapBytes4IfLittleEndian (p);
932     double dMaxValue = *reinterpret_cast<float*>(p);
933
934     DetectorArray& detArray = getDetectorArray (iv);
935     detArray.setViewAngle (dAlpha);
936     DetectorValue* detval = detArray.detValues();
937
938     double dViewScale = 1. / (2294.4871 * ::pow (2.0, -lEscale));
939     lDataPos += 32;
940     for (int id = 0; id < iNDets; id++) {
941       int iV = pData[lDataPos+1] + (pData[lDataPos] << 8);
942       if (iV > 32767)   // two's complement signed conversion
943         iV = iV - 65536;
944       detval[id] = iV * dViewScale * pdCosScale[id];
945       lDataPos += 2;
946     }
947   }
948
949   delete pdCosScale;
950   return true;
951 }
952
953 Projections*
954 Projections::interpolateToParallel () const
955 {
956   if (m_geometry == Scanner::GEOMETRY_PARALLEL)
957     return const_cast<Projections*>(this);
958
959   int nDet = m_nDet;
960   int nView = m_nView;
961   Projections* pProjNew = new Projections (nView, nDet);
962   pProjNew->m_geometry = Scanner::GEOMETRY_PARALLEL;
963   pProjNew->m_dFocalLength = m_dFocalLength;
964   pProjNew->m_dSourceDetectorLength = m_dSourceDetectorLength;
965   pProjNew->m_dViewDiameter = m_dViewDiameter;
966   pProjNew->m_dFanBeamAngle = m_dFanBeamAngle;
967   pProjNew->m_calcTime  = 0;
968   pProjNew->m_remark = m_remark;
969   pProjNew->m_remark += "; Interpolate to Parallel";
970   pProjNew->m_label.setLabelType (Array2dFileLabel::L_HISTORY);
971   pProjNew->m_label.setLabelString (pProjNew->m_remark);
972   pProjNew->m_label.setCalcTime (pProjNew->m_calcTime);
973   pProjNew->m_label.setDateTime (pProjNew->m_year, pProjNew->m_month, pProjNew->m_day, pProjNew->m_hour, pProjNew->m_minute, pProjNew->m_second);
974
975   pProjNew->m_rotStart = 0;
976 #ifdef CONVERT_PARALLEL_PI
977   pProjNew->m_rotInc = PI / nView;;
978 #else
979   pProjNew->m_rotInc = TWOPI / nView;
980 #endif
981   pProjNew->m_detStart = -m_dViewDiameter / 2;
982   pProjNew->m_detInc = m_dViewDiameter / nDet;
983   if (nDet % 2 == 0) // even
984     pProjNew->m_detInc = m_dViewDiameter / (nDet - 1);
985
986   ParallelRaysums parallel (this, ParallelRaysums::THETA_RANGE_NORMALIZE_TO_TWOPI);
987
988   double* pdThetaValuesForT = new double [pProjNew->nView()];
989   double* pdRaysumsForT = new double [pProjNew->nView()];
990
991   // interpolate to evenly spaced theta (views)
992   double dDetPos = pProjNew->m_detStart;
993   for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc) {
994     parallel.getThetaAndRaysumsForT (iD, pdThetaValuesForT, pdRaysumsForT);
995
996     double dViewAngle = m_rotStart;
997     int iLastFloor = -1;
998     for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
999       DetectorValue* detValues = pProjNew->getDetectorArray (iV).detValues();
1000
1001       detValues[iD] = parallel.interpolate (pdThetaValuesForT, pdRaysumsForT, pProjNew->nView(), dViewAngle, &iLastFloor);
1002     }
1003   }
1004   delete pdThetaValuesForT;
1005   delete pdRaysumsForT;
1006
1007   // interpolate to evenly space t (detectors)
1008   double* pdOriginalDetPositions = new double [pProjNew->nDet()];
1009   parallel.getDetPositions (pdOriginalDetPositions);
1010
1011   double* pdDetValueCopy = new double [pProjNew->nDet()];
1012   double dViewAngle = m_rotStart;
1013   for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
1014     DetectorArray& detArray = pProjNew->getDetectorArray (iV);
1015     DetectorValue* detValues = detArray.detValues();
1016     detArray.setViewAngle (dViewAngle);
1017
1018     for (int i = 0; i < pProjNew->nDet(); i++)
1019       pdDetValueCopy[i] = detValues[i];
1020
1021     double dDetPos = pProjNew->m_detStart;
1022     int iLastFloor = -1;
1023     for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc) {
1024       detValues[iD] = parallel.interpolate (pdOriginalDetPositions, pdDetValueCopy, pProjNew->nDet(), dDetPos, &iLastFloor);
1025     }
1026   }
1027   delete pdDetValueCopy;
1028   delete pdOriginalDetPositions;
1029
1030   return pProjNew;
1031 }
1032
1033
1034 ///////////////////////////////////////////////////////////////////////////////
1035 //
1036 // Class ParallelRaysums
1037 //
1038 // Used for converting divergent beam raysums into Parallel raysums
1039 //
1040 ///////////////////////////////////////////////////////////////////////////////
1041
1042 ParallelRaysums::ParallelRaysums (const Projections* pProjections, int iThetaRange)
1043 : m_iNumCoordinates(0), m_iNumView(pProjections->nView()), m_iNumDet(pProjections->nDet()),
1044   m_iThetaRange (iThetaRange), m_pCoordinates(NULL)
1045 {
1046   int iGeometry = pProjections->geometry();
1047   double dDetInc = pProjections->detInc();
1048   double dDetStart = pProjections->detStart();
1049   double dFocalLength = pProjections->focalLength();
1050
1051   m_iNumCoordinates =  m_iNumView * m_iNumDet;
1052   m_pCoordinates = new ParallelRaysumCoordinate [m_iNumCoordinates];
1053   m_vecpCoordinates.reserve (m_iNumCoordinates);
1054   for (int i = 0; i < m_iNumCoordinates; i++)
1055     m_vecpCoordinates[i] = m_pCoordinates + i;
1056
1057   int iCoordinate = 0;
1058   for (int iV = 0; iV < m_iNumView; iV++) {
1059     double dViewAngle = pProjections->getDetectorArray(iV).viewAngle();
1060     const DetectorValue* detValues = pProjections->getDetectorArray(iV).detValues();
1061
1062     double dDetPos = dDetStart;
1063     for (int iD = 0; iD < m_iNumDet; iD++) {
1064       ParallelRaysumCoordinate* pC = m_vecpCoordinates[iCoordinate++];
1065
1066       if (iGeometry == Scanner::GEOMETRY_PARALLEL) {
1067         pC->m_dTheta = dViewAngle;
1068         pC->m_dT = dDetPos;
1069       } else if (iGeometry == Scanner::GEOMETRY_EQUILINEAR) {
1070         double dFanAngle = atan (dDetPos / pProjections->sourceDetectorLength());
1071         pC->m_dTheta = dViewAngle + dFanAngle;
1072         pC->m_dT = dFocalLength * sin(dFanAngle);        
1073
1074       } else if (iGeometry == Scanner::GEOMETRY_EQUIANGULAR) {
1075         // fan angle is same as dDetPos
1076         pC->m_dTheta = dViewAngle + dDetPos;
1077         pC->m_dT = dFocalLength * sin (dDetPos);        
1078       }
1079       if (m_iThetaRange != THETA_RANGE_UNCONSTRAINED) {
1080         pC->m_dTheta = normalizeAngle (pC->m_dTheta);
1081         if (m_iThetaRange == THETA_RANGE_FOLD_TO_PI && pC->m_dTheta >= PI) {
1082           pC->m_dTheta -= PI;
1083           pC->m_dT = -pC->m_dT;
1084         }
1085       }
1086       pC->m_dRaysum = detValues[iD];
1087       dDetPos += dDetInc;
1088     }
1089   }
1090 }
1091
1092 ParallelRaysums::~ParallelRaysums()
1093 {
1094   delete m_pCoordinates;
1095 }
1096
1097 ParallelRaysums::CoordinateContainer&
1098 ParallelRaysums::getSortedByTheta()
1099 {
1100   if (m_vecpSortedByTheta.size() == 0) {
1101     m_vecpSortedByTheta.resize (m_iNumCoordinates);
1102     for (int i = 0; i < m_iNumCoordinates; i++)
1103       m_vecpSortedByTheta[i] = m_vecpCoordinates[i];
1104     std::sort (m_vecpSortedByTheta.begin(), m_vecpSortedByTheta.end(), ParallelRaysumCoordinate::compareByTheta);
1105   }
1106
1107   return m_vecpSortedByTheta;
1108 }
1109
1110 ParallelRaysums::CoordinateContainer&
1111 ParallelRaysums::getSortedByT()
1112 {
1113   if (m_vecpSortedByT.size() == 0) {
1114     m_vecpSortedByT.resize (m_iNumCoordinates);
1115     for (int i = 0; i < m_iNumCoordinates; i++)
1116       m_vecpSortedByT[i] = m_vecpCoordinates[i];
1117     std::sort (m_vecpSortedByT.begin(), m_vecpSortedByT.end(), ParallelRaysumCoordinate::compareByT);
1118   }
1119
1120   return m_vecpSortedByT;
1121 }
1122
1123
1124 void
1125 ParallelRaysums::getLimits (double* dMinT, double* dMaxT, double* dMinTheta, double* dMaxTheta) const
1126 {
1127   if (m_iNumCoordinates <= 0)
1128     return;
1129
1130   *dMinT = *dMaxT = m_vecpCoordinates[0]->m_dT;
1131   *dMinTheta = *dMaxTheta = m_vecpCoordinates[0]->m_dTheta;
1132
1133   for (int i = 0; i < m_iNumCoordinates; i++) {
1134     double dT = m_vecpCoordinates[i]->m_dT;
1135     double dTheta = m_vecpCoordinates[i]->m_dTheta;
1136
1137     if (dT < *dMinT)
1138       *dMinT = dT;
1139     else if (dT > *dMaxT)
1140       *dMaxT = dT;
1141
1142     if (dTheta < *dMinTheta)
1143       *dMinTheta = dTheta;
1144     else if (dTheta > *dMaxTheta)
1145       *dMaxTheta = dTheta;
1146   }
1147 }
1148
1149 void
1150 ParallelRaysums::getThetaAndRaysumsForT (int iTheta, double* pTheta, double* pRaysum)
1151 {
1152   const CoordinateContainer& coordsT = getSortedByT();
1153
1154   int iBase = iTheta * m_iNumView;
1155   for (int i = 0; i < m_iNumView; i++) {
1156     int iPos = iBase + i;
1157     pTheta[i] = coordsT[iPos]->m_dTheta;
1158     pRaysum[i] = coordsT[iPos]->m_dRaysum;
1159   }
1160 }
1161
1162 void
1163 ParallelRaysums::getDetPositions (double* pdDetPos)
1164 {
1165   const CoordinateContainer& coordsT = getSortedByT();
1166
1167   int iPos = 0;
1168   for (int i = 0; i < m_iNumDet; i++) {
1169     pdDetPos[i] = coordsT[iPos]->m_dT;
1170     iPos += m_iNumView;
1171   }
1172 }
1173
1174 // locate by bisection, O(log2(n))
1175 // iLastFloor is used when sequential calls to interpolate have monotonically increasing dX
1176 double
1177 ParallelRaysums::interpolate (double* pdX, double* pdY, int n, double dX, int* iLastFloor)
1178 {
1179   int iLower = -1;
1180   int iUpper = n;
1181   if (iLastFloor && *iLastFloor >= 0 && pdX[*iLastFloor] < dX)
1182     iLower = *iLastFloor;
1183
1184   while (iUpper - iLower > 1) {
1185     int iMiddle = (iUpper + iLower) >> 1;
1186     if (dX >= pdX[iMiddle])
1187       iLower = iMiddle;
1188     else
1189       iUpper = iMiddle;
1190   }
1191   if (dX <= pdX[0])
1192     return pdY[0];
1193   else if (dX >= pdX[n-1])
1194     return pdY[1];
1195
1196   if (iLower < 0 || iLower >= n) {
1197     sys_error (ERR_SEVERE, "Coordinate out of range [locateThetaBase]");
1198     return 0;
1199   }
1200
1201   if (iLastFloor)
1202     *iLastFloor = iLower;
1203   return pdY[iLower] + (pdY[iUpper] - pdY[iLower]) * ((dX - pdX[iLower]) / (pdX[iUpper] - pdX[iLower]));
1204 }
1205