d80579e5f897ab0ea3d17372dd5ed6d5fd5fecdd
[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.48 2001/02/20 17:44:14 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* Projections::s_aszInterpName[] = 
37 {
38   {"nearest"},
39   {"bilinear"},
40 //  {"bicubic"},
41 };
42
43 const char* 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_geometry = scanner.geometry();
154   m_dFocalLength = scanner.focalLength();
155   m_dViewDiameter = scanner.viewDiameter();
156   m_rotStart = 0;
157   m_detStart =  -(scanner.detLen() / 2);
158   m_dFanBeamAngle = scanner.fanBeamAngle();
159 }
160
161 void
162 Projections::setNView (int nView)  // used by MPI to reduce # of views
163 {
164   deleteProjData();
165   init (nView, m_nDet);
166 }
167
168 // NAME
169 // newProjData
170
171 void 
172 Projections::newProjData (void)
173 {
174   if (m_projData)
175     sys_error(ERR_WARNING, "m_projData != NULL [newProjData]");
176   
177   if (m_nView > 0 && m_nDet) {
178     m_projData = new DetectorArray* [m_nView];
179     
180     for (int i = 0; i < m_nView; i++)
181       m_projData[i] = new DetectorArray (m_nDet);
182   }
183 }
184
185
186 /* NAME
187 *   projections_free                    Free memory allocated to projections
188 *
189 * SYNOPSIS
190 *   projections_free(proj)
191 *   Projections& proj                           Projectionss to be deallocated
192 */
193
194 void 
195 Projections::deleteProjData (void)
196 {
197   if (m_projData != NULL) {
198     for (int i = 0; i < m_nView; i++)
199       delete m_projData[i];
200     
201     delete m_projData;
202     m_projData = NULL;
203   }
204 }
205
206
207 /* NAME
208 *    Projections::headerWwrite         Write data header for projections file
209 *
210 */
211
212 bool 
213 Projections::headerWrite (fnetorderstream& fs)
214 {
215   kuint16 _hsize = m_headerSize;
216   kuint16 _signature = m_signature;
217   kuint32 _nView = m_nView;
218   kuint32 _nDet = m_nDet;
219   kuint32 _geom = m_geometry;
220   kuint16 _remarksize = m_remark.length();
221   kuint16 _year = m_year;
222   kuint16 _month = m_month;
223   kuint16 _day = m_day;
224   kuint16 _hour = m_hour;
225   kuint16 _minute = m_minute;
226   kuint16 _second = m_second;
227   
228   kfloat64 _calcTime = m_calcTime;
229   kfloat64 _rotStart = m_rotStart;
230   kfloat64 _rotInc = m_rotInc;
231   kfloat64 _detStart = m_detStart;
232   kfloat64 _detInc = m_detInc;
233   kfloat64 _viewDiameter = m_dViewDiameter;
234   kfloat64 _focalLength = m_dFocalLength;
235   kfloat64 _fanBeamAngle = m_dFanBeamAngle;
236
237   fs.seekp(0);
238   if (! fs)
239     return false;
240   
241   fs.writeInt16 (_hsize);
242   fs.writeInt16 (_signature);
243   fs.writeInt32 (_nView);
244   fs.writeInt32 (_nDet);
245   fs.writeInt32 (_geom);
246   fs.writeFloat64 (_calcTime);
247   fs.writeFloat64 (_rotStart);
248   fs.writeFloat64 (_rotInc);
249   fs.writeFloat64 (_detStart);
250   fs.writeFloat64 (_detInc);
251   fs.writeFloat64 (_viewDiameter);
252   fs.writeFloat64 (_focalLength);
253   fs.writeFloat64 (_fanBeamAngle);
254   fs.writeInt16 (_year);
255   fs.writeInt16 (_month);
256   fs.writeInt16 (_day);
257   fs.writeInt16 (_hour);
258   fs.writeInt16 (_minute);
259   fs.writeInt16 (_second);
260   fs.writeInt16 (_remarksize);
261   fs.write (m_remark.c_str(), _remarksize);
262   
263   m_headerSize = fs.tellp();
264   _hsize = m_headerSize;
265   fs.seekp(0);
266   fs.writeInt16 (_hsize);
267   if (! fs)
268     return false;
269   
270   return true;
271 }
272
273 /* NAME
274 *    projections_read_header         Read data header for projections file
275 *
276 */
277 bool
278 Projections::headerRead (fnetorderstream& fs)
279 {
280   kuint16 _hsize, _signature, _year, _month, _day, _hour, _minute, _second, _remarksize = 0;
281   kuint32 _nView, _nDet, _geom;
282   kfloat64 _calcTime, _rotStart, _rotInc, _detStart, _detInc, _focalLength, _viewDiameter, _fanBeamAngle;
283   
284   fs.seekg(0);
285   if (! fs)
286     return false;
287   
288   fs.readInt16 (_hsize);
289   fs.readInt16 (_signature);
290   fs.readInt32 (_nView);
291   fs.readInt32 (_nDet);
292   fs.readInt32 (_geom);
293   fs.readFloat64 (_calcTime);
294   fs.readFloat64 (_rotStart);
295   fs.readFloat64 (_rotInc);
296   fs.readFloat64 (_detStart);
297   fs.readFloat64 (_detInc);
298   fs.readFloat64 (_viewDiameter);
299   fs.readFloat64 (_focalLength);
300   fs.readFloat64 (_fanBeamAngle);
301   fs.readInt16 (_year);
302   fs.readInt16 (_month);
303   fs.readInt16 (_day);
304   fs.readInt16 (_hour);
305   fs.readInt16 (_minute);
306   fs.readInt16 (_second);
307   fs.readInt16 (_remarksize);
308   
309   if (! fs) {
310     sys_error (ERR_SEVERE, "Error reading header information , _remarksize=%d [projections_read_header]", _remarksize);
311     return false;
312   }
313   
314   if (_signature != m_signature) {
315     sys_error (ERR_SEVERE, "File %s does not have a valid projection file signature", m_filename.c_str());
316     return false;
317   }
318   
319   char* pszRemarkStorage = new char [_remarksize+1];
320   fs.read (pszRemarkStorage, _remarksize);
321   if (! fs) {
322     sys_error (ERR_SEVERE, "Error reading remark, _remarksize = %d", _remarksize);
323     return false;
324   }
325   pszRemarkStorage[_remarksize] = 0;
326   m_remark = pszRemarkStorage;
327   delete pszRemarkStorage;
328   
329   off_t _hsizeread = fs.tellg();
330   if (!fs || _hsizeread != _hsize) {
331     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);
332     return false;
333   }
334   
335   m_headerSize = _hsize;
336   m_nView = _nView;
337   m_nDet = _nDet;
338   m_geometry = _geom;
339   m_calcTime = _calcTime;
340   m_rotStart = _rotStart;
341   m_rotInc = _rotInc;
342   m_detStart = _detStart;
343   m_detInc = _detInc;
344   m_dFocalLength = _focalLength;
345   m_dViewDiameter = _viewDiameter;
346   m_dFanBeamAngle = _fanBeamAngle;
347   m_year = _year;
348   m_month = _month;
349   m_day = _day;
350   m_hour = _hour;
351   m_minute = _minute;
352   m_second = _second;
353   
354   m_label.setLabelType (Array2dFileLabel::L_HISTORY);
355   m_label.setLabelString (m_remark);
356   m_label.setCalcTime (m_calcTime);
357   m_label.setDateTime (m_year, m_month, m_day, m_hour, m_minute, m_second);
358   
359   return true;
360 }
361
362 bool
363 Projections::read (const std::string& filename)
364 {
365   return read (filename.c_str());
366 }
367
368
369 bool
370 Projections::read (const char* filename)
371 {
372   m_filename = filename;
373 #ifdef MSVC
374   frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary);
375 #else
376   frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary | std::ios::nocreate);
377 #endif
378   
379   if (fileRead.fail())
380     return false;
381   
382   if (! headerRead (fileRead))
383     return false;
384   
385   deleteProjData ();
386   newProjData();
387   
388   for (int i = 0; i < m_nView; i++) {
389     if (! detarrayRead (fileRead, *m_projData[i], i))
390       break;
391   }
392   
393   fileRead.close();
394   return true;
395 }
396
397
398 bool 
399 Projections::copyViewData (const std::string& filename, std::ostream& os, int startView, int endView)
400 {
401   return copyViewData (filename.c_str(), os, startView, endView);
402 }
403
404 bool 
405 Projections::copyViewData (const char* const filename, std::ostream& os, int startView, int endView)
406 {
407   frnetorderstream is (filename, std::ios::in | std::ios::binary);
408   kuint16 sizeHeader, signature;
409   kuint32 _nView, _nDet;
410   
411   is.seekg (0);
412   if (is.fail()) {
413     sys_error (ERR_SEVERE, "Unable to read projection file %s", filename);
414     return false;
415   }
416
417   is.readInt16 (sizeHeader);
418   is.readInt16 (signature);
419   is.readInt32 (_nView);
420   is.readInt32 (_nDet);
421   int nView = _nView;
422   int nDet = _nDet;
423   
424   if (signature != m_signature) {
425     sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
426     return false;
427   }
428   
429   if (startView < 0)
430     startView = 0;
431   if (startView > nView - 1)
432     startView = nView;
433   if (endView < 0 || endView > nView - 1)
434     endView = nView - 1;
435   
436   if (startView > endView) { // swap if start > end
437     int tempView = endView;
438     endView = startView;
439     startView = tempView;
440   }
441   
442   int sizeView = 8 /* view_angle */ + 4 /* nDet */ + (4 * nDet);
443   unsigned char* pViewData = new unsigned char [sizeView];
444   
445   for (int i = startView; i <= endView; i++) {
446     is.seekg (sizeHeader + i * sizeView);
447     is.read (reinterpret_cast<char*>(pViewData), sizeView);
448     os.write (reinterpret_cast<char*>(pViewData), sizeView);
449     if (is.fail() || os.fail())
450       break;
451   }
452   
453   delete pViewData;
454   if (is.fail()) 
455     sys_error (ERR_SEVERE, "Error reading projection file");
456   if (os.fail()) 
457     sys_error (ERR_SEVERE, "Error writing projection file");
458   
459   return (! (is.fail() | os.fail()));
460 }
461
462 bool 
463 Projections::copyHeader (const std::string& filename, std::ostream& os)
464 {
465   return copyHeader (filename.c_str(), os);
466 }
467
468 bool
469 Projections::copyHeader (const char* const filename, std::ostream& os)
470 {
471   frnetorderstream is (filename, std::ios::in | std::ios::binary);
472   kuint16 sizeHeader, signature;
473   is.readInt16 (sizeHeader);
474   is.readInt16 (signature);
475   is.seekg (0);
476   if (signature != m_signature) {
477     sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
478     return false;
479   }
480   
481   unsigned char* pHdrData = new unsigned char [sizeHeader];
482   is.read (reinterpret_cast<char*>(pHdrData), sizeHeader);
483   if (is.fail()) {
484     sys_error (ERR_SEVERE, "Error reading header");
485     return false;
486   }
487   
488   os.write (reinterpret_cast<char*>(pHdrData), sizeHeader);
489   if (os.fail()) {
490     sys_error (ERR_SEVERE, "Error writing header");
491     return false;
492   }
493   
494   return true;
495 }
496
497 bool
498 Projections::write (const std::string& filename)
499 {
500   return write (filename.c_str());
501 }
502
503 bool
504 Projections::write (const char* filename)
505 {
506   frnetorderstream fs (filename, std::ios::out | std::ios::binary | std::ios::trunc | std::ios::ate);
507   m_filename = filename;
508   if (! fs) {
509     sys_error (ERR_SEVERE, "Error opening file %s for output [projections_create]", filename);
510     return false;
511   }
512   
513   if (! headerWrite (fs))
514     return false;
515   
516   if (m_projData != NULL) {
517     for (int i = 0; i < m_nView; i++) {
518       if (! detarrayWrite (fs, *m_projData[i], i))
519         break;
520     }
521   }
522   if (! fs)
523     return false;
524   
525   fs.close();
526   
527   return true;
528 }
529
530 /* NAME
531 *   detarrayRead                Read a Detector Array structure from the disk
532 *
533 * SYNOPSIS
534 *   detarrayRead (proj, darray, view_num)
535 *   DETARRAY *darray            Detector array storage location to be filled
536 *   int      view_num           View number to read
537 */
538
539 bool
540 Projections::detarrayRead (fnetorderstream& fs, DetectorArray& darray, const int iview)
541 {
542   const int detval_bytes = darray.nDet() * sizeof(kfloat32);
543   const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
544   const int view_bytes = detheader_bytes + detval_bytes;
545   const off_t start_data = m_headerSize + (iview * view_bytes);
546   DetectorValue* detval_ptr = darray.detValues();  
547   kfloat64 view_angle;
548   kuint32 nDet;
549   
550   fs.seekg (start_data);
551   
552   fs.readFloat64 (view_angle);
553   fs.readInt32 (nDet);
554   darray.setViewAngle (view_angle);
555   //  darray.setNDet ( nDet);
556   
557   for (unsigned int i = 0; i < nDet; i++) {
558     kfloat32 detval;
559     fs.readFloat32 (detval);
560     detval_ptr[i] = detval;
561   }
562   if (! fs)
563     return false;
564   
565   return true;
566 }
567
568
569 /* NAME
570 *   detarrayWrite                       Write detector array data to the disk
571 *
572 * SYNOPSIS
573 *   detarrayWrite (darray, view_num)
574 *   DETARRAY *darray                    Detector array data to be written
575 *   int      view_num                   View number to write
576 *
577 * DESCRIPTION
578 *       This routine writes the detarray data from the disk sequentially to
579 *    the file that was opened with open_projections().  Data is written in
580 *    binary format.
581 */
582
583 bool
584 Projections::detarrayWrite (fnetorderstream& fs, const DetectorArray& darray, const int iview)
585 {
586   const int detval_bytes = darray.nDet() * sizeof(float);
587   const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
588   const int view_bytes = detheader_bytes + detval_bytes;
589   const off_t start_data = m_headerSize + (iview * view_bytes);
590   const DetectorValue* const detval_ptr = darray.detValues();  
591   kfloat64 view_angle = darray.viewAngle();
592   kuint32 nDet = darray.nDet();
593   
594   fs.seekp (start_data);
595   if (! fs) {
596     sys_error (ERR_SEVERE, "Error seeking detectory array [detarrayWrite]");
597     return false;
598   }
599   
600   fs.writeFloat64 (view_angle);
601   fs.writeInt32 (nDet);
602   
603   for (unsigned int i = 0; i < nDet; i++) {
604     kfloat32 detval = detval_ptr[i];
605     fs.writeFloat32 (detval);
606   }
607   
608   if (! fs)
609     return (false);
610   
611   return true;
612 }
613
614 /* NAME
615 *   printProjectionData                 Print projections data
616 *
617 * SYNOPSIS
618 *   printProjectionData ()
619 */
620
621 void
622 Projections::printProjectionData ()
623 {
624   printProjectionData (0, nView() - 1);
625 }
626
627 void
628 Projections::printProjectionData (int startView, int endView)
629 {
630   printf("Projections Data\n\n");
631   printf("Description: %s\n", m_remark.c_str());
632   printf("Geometry: %s\n", Scanner::convertGeometryIDToName (m_geometry));
633   printf("nView       = %8d           nDet = %8d\n", m_nView, m_nDet);
634   printf("focalLength = %8.4f ViewDiameter = %8.4f\n", m_dFocalLength, m_dViewDiameter);
635   printf("fanBeamAngle= %8.4f\n", convertRadiansToDegrees(m_dFanBeamAngle));
636   printf("rotStart    = %8.4f       rotInc = %8.4f\n", m_rotStart, m_rotInc);
637   printf("detStart    = %8.4f       detInc = %8.4f\n", m_detStart, m_detInc);
638   if (m_projData != NULL) {
639     if (startView < 0)
640       startView = 0;
641     if (endView < 0)
642       endView = m_nView - 1;
643     if (startView > m_nView - 1)
644       startView = m_nView - 1;
645     if (endView > m_nView - 1)
646       endView = m_nView - 1;
647     for (int ir = startView; ir <= endView - 1; ir++) {
648       printf("View %d: angle %f\n", ir, m_projData[ir]->viewAngle());
649       DetectorValue* detval = m_projData[ir]->detValues();
650       for (int id = 0; id < m_projData[ir]->nDet(); id++)
651         printf("%8.4f  ", detval[id]);
652       printf("\n");
653     }
654   }
655 }
656
657 void 
658 Projections::printScanInfo (std::ostringstream& os) const
659 {
660   os << "Number of detectors: " << m_nDet << "\n";
661   os << "Number of views: " << m_nView<< "\n";
662   os << "Description: " << m_remark.c_str()<< "\n";
663   os << "Geometry: " << Scanner::convertGeometryIDToName (m_geometry)<< "\n";
664   os << "Focal Length: " << m_dFocalLength<< "\n";
665   os << "View Diameter: " << m_dViewDiameter<< "\n";
666   os << "Fan Beam Angle: " << convertRadiansToDegrees(m_dFanBeamAngle) << "\n";
667   os << "detStart: " << m_detStart<< "\n";
668   os << "detInc: " << m_detInc<< "\n";
669   os << "rotStart: " << m_rotStart<< "\n";
670   os << "rotInc: " << m_rotInc<< "\n";
671 }
672
673
674 bool 
675 Projections::convertPolar (ImageFile& rIF, int iInterpolationID)
676 {
677   unsigned int nx = rIF.nx();
678   unsigned int ny = rIF.ny();
679   ImageFileArray v = rIF.getArray();
680   ImageFileArray vImag = rIF.getImaginaryArray();
681
682   if (! v || nx == 0 || ny == 0)
683     return false;
684   
685   Array2d<double> adView (nx, ny);
686   Array2d<double> adDet (nx, ny);
687   double** ppdView = adView.getArray();
688   double** ppdDet = adDet.getArray();
689
690   if (! calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet)) 
691     return false;
692
693   std::complex<double>** ppcDetValue = new std::complex<double>* [m_nView];
694   unsigned int iView;
695   for (iView = 0; iView < m_nView; iView++) {
696     ppcDetValue[iView] = new std::complex<double> [m_nDet];
697     for (unsigned int iDet = 0; iDet < m_nDet; iDet++)
698       ppcDetValue[iView][iDet] = std::complex<double>(getDetectorArray (iView).detValues()[iDet], 0);
699   }
700
701   interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, m_nView, m_nDet, iInterpolationID);
702
703   for (iView = 0; iView < m_nView; iView++)
704     delete [] ppcDetValue[iView];
705   delete [] ppcDetValue;
706
707   return true;
708 }
709
710
711 bool 
712 Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad)
713 {
714   unsigned int nx = rIF.nx();
715   unsigned int ny = rIF.ny();
716   ImageFileArray v = rIF.getArray();
717   if (! rIF.isComplex())
718     rIF.convertRealToComplex();
719   ImageFileArray vImag = rIF.getImaginaryArray();
720
721   if (! v || nx == 0 || ny == 0)
722     return false;
723   
724 #ifndef HAVE_FFT
725   return false;
726 #else
727   Array2d<double> adView (nx, ny);
728   Array2d<double> adDet (nx, ny);
729   double** ppdView = adView.getArray();
730   double** ppdDet = adDet.getArray();
731
732   std::complex<double>** ppcDetValue = new std::complex<double>* [m_nView];
733   unsigned int iView;
734   double* pdDet = new double [m_nDet];
735   fftw_complex* pcIn = new fftw_complex [m_nDet];
736   fftw_plan plan = fftw_create_plan (m_nDet, FFTW_FORWARD, FFTW_IN_PLACE);
737
738   for (iView = 0; iView < m_nView; iView++) {
739     unsigned int iDet;
740     for (iDet = 0; iDet < m_nDet; iDet++) {
741       pcIn[iDet].re = getDetectorArray(iView).detValues()[iDet];
742       pcIn[iDet].im = 0;
743     }
744     fftw_one (plan, pcIn, NULL);
745     ppcDetValue[iView] = new std::complex<double> [m_nDet];
746     for (iDet = 0; iDet < m_nDet; iDet++)
747       ppcDetValue[iView][iDet] = std::complex<double> (pcIn[iDet].re, pcIn[iDet].im); 
748     Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], m_nDet);
749   }
750
751   fftw_destroy_plan (plan);  
752   delete [] pcIn;
753   
754   bool bError = calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet);
755
756   if (! bError)
757     interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, m_nView, m_nDet, iInterpolationID);
758
759   for (iView = 0; iView < m_nView; iView++)
760     delete [] ppcDetValue[iView];
761   delete [] ppcDetValue;
762
763   return bError;
764 #endif
765 }
766
767
768 bool
769 Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double** ppdView, double** ppdDet)
770 {
771   double xMin = -phmLen() / 2;
772   double xMax = xMin + phmLen();
773   double yMin = -phmLen() / 2;
774   double yMax = yMin + phmLen();
775   
776   double xInc = (xMax - xMin) / nx;     // size of cells
777   double yInc = (yMax - yMin) / ny;
778   
779   int iDetCenter = (m_nDet - 1) / 2;    // index refering to L=0 projection 
780
781   if (m_geometry != Scanner::GEOMETRY_PARALLEL) {
782     sys_error (ERR_WARNING, "convertPolar supports Parallel only");
783     return false;
784   }
785   
786   // Calculates polar coordinates (view#, det#) for each point on phantom grid
787   double x = xMin + xInc / 2;   // Rectang coords of center of pixel 
788   for (unsigned int ix = 0; ix < nx; x += xInc, ix++) {
789     double y = yMin + yInc / 2;
790     for (unsigned int iy = 0; iy < ny; y += yInc, iy++) {
791       double r = ::sqrt (x * x + y * y);
792       double phi = atan2 (y, x);
793
794       if (phi >= PI) {
795         phi -= PI;
796       } else if (phi < 0) {
797         phi += PI;
798       } else
799         r = -r;
800       
801       ppdView[ix][iy] = (phi - m_rotStart) / m_rotInc;
802       ppdDet[ix][iy] = (r / m_detInc) + iDetCenter;
803     }
804   }
805
806   return true;
807 }
808
809 void
810 Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag,
811      unsigned int nx, unsigned int ny, std::complex<double>** ppcDetValue,
812      double** ppdView, double** ppdDet, unsigned int nView, unsigned int nDet, int iInterpolationID)
813 {
814   for (unsigned int ix = 0; ix < ny; ix++) {
815     for (unsigned int iy = 0; iy < ny; iy++) {
816       if (iInterpolationID == POLAR_INTERP_NEAREST) {
817         int iView = nearest<int> (ppdView[ix][iy]);
818         int iDet = nearest<int> (ppdDet[ix][iy]);
819         if (iView == nView) {
820           iView = 0;
821        //   iDet = m_nDet - iDet;
822         }
823         if (iDet >= 0 && iDet < nDet && iView >= 0 && iView < nView) {
824           v[ix][iy] = ppcDetValue[iView][iDet].real();
825           if (vImag)
826             vImag[ix][iy] = ppcDetValue[iView][iDet].imag();
827         } else {
828           sys_error (ERR_SEVERE, "Can't find projection data for ix=%d,iy=%d with radView=%f and radDet=%f", 
829             ix, iy, ppdView[ix][iy], ppdDet[ix][iy]);
830           v[ix][iy] = 0;
831         }
832       } else if (iInterpolationID == POLAR_INTERP_BILINEAR) {
833         int iFloorView = static_cast<int>(ppdView[ix][iy]);
834         double dFracView = ppdView[ix][iy] - iFloorView;
835         int iFloorDet = static_cast<int>(ppdDet[ix][iy]);
836         double dFracDet = ppdDet[ix][iy] - iFloorDet;
837
838         if (iFloorDet >= 0 && iFloorView >= 0) { 
839           std::complex<double> v1 = ppcDetValue[iFloorView][iFloorDet];
840           std::complex<double> v2, v3, v4;
841           if (iFloorView < nView - 1)
842             v2 = ppcDetValue[iFloorView + 1][iFloorDet];
843           else 
844             v2 = ppcDetValue[0][iFloorDet];
845           if (iFloorDet < nDet - 1) 
846             v4 = ppcDetValue[iFloorView][iFloorDet+1];
847           else
848             v4 = v1;
849           if (iFloorView < nView - 1 && iFloorDet < nDet - 1)
850             v3 = ppcDetValue [iFloorView+1][iFloorDet+1];
851           else if (iFloorView < nView - 1)
852             v3 = v2;
853           else 
854             v3 = ppcDetValue[0][iFloorDet+1];
855           std::complex<double> vInterp = (1 - dFracView) * (1 - dFracDet) * v1 +
856             dFracView * (1 - dFracDet) * v2 + dFracView * dFracDet * v3 +
857             dFracDet * (1 - dFracView) * v4;
858           v[ix][iy] = vInterp.real();
859           if (vImag)
860             vImag[ix][iy] = vInterp.imag();
861         } else {
862           sys_error (ERR_SEVERE, "Can't find projection data for ix=%d,iy=%d with radView=%f and radDet=%f", 
863             ix, iy, ppdView[ix][iy], ppdDet[ix][iy]);
864           v[ix][iy] = 0;
865           if (vImag)
866             vImag[ix][iy] = 0;
867         }
868       } else if (iInterpolationID == POLAR_INTERP_BICUBIC) {
869         v[ix][iy] =0;
870           if (vImag)
871             vImag[ix][iy] = 0;
872       }
873     }
874   }
875 }
876
877