r355: Polar conversions of projections
[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-2000 Kevin Rosenberg
10 **
11 **  $Id: projections.cpp,v 1.38 2001/01/06 15:33:15 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_phmLen = scanner.phmLen();
152   m_rotInc = scanner.rotInc();
153   m_detInc = scanner.detInc();
154   m_geometry = scanner.geometry();
155   m_focalLength = scanner.focalLength();
156   m_fieldOfView = scanner.fieldOfView();
157   m_rotStart = 0;
158   m_detStart =  -(scanner.detLen() / 2);
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 _phmLen = m_phmLen;
234   kfloat64 _fieldOfView = m_fieldOfView;
235   kfloat64 _focalLength = m_focalLength;
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 (_phmLen);
252   fs.writeFloat64 (_focalLength);
253   fs.writeFloat64 (_fieldOfView);
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, _phmLen, _focalLength, _fieldOfView;
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 (_phmLen);
299   fs.readFloat64 (_focalLength);
300   fs.readFloat64 (_fieldOfView);
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_phmLen = _phmLen;
345   m_focalLength = _focalLength;
346   m_fieldOfView = _fieldOfView;
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, ios::in | ios::binary);
408   kuint16 sizeHeader, signature;
409   kuint32 _nView, _nDet;
410   
411   is.readInt16 (sizeHeader);
412   is.readInt16 (signature);
413   is.readInt32 (_nView);
414   is.readInt32 (_nDet);
415   int nView = _nView;
416   int nDet = _nDet;
417   
418   if (signature != m_signature) {
419     sys_error (ERR_FATAL, "Illegal signature in projection file %s", filename);
420     return false;
421   }
422   
423   if (startView < 0)
424     startView = 0;
425   if (startView > nView - 1)
426     startView = nView;
427   if (endView < 0 || endView > nView - 1)
428     endView = nView - 1;
429   
430   if (startView > endView) { // swap if start > end
431     int tempView = endView;
432     endView = startView;
433     startView = tempView;
434   }
435   
436   int sizeView = 8 /* view_angle */ + 4 /* nDet */ + (4 * nDet);
437   unsigned char* pViewData = new unsigned char [sizeView];
438   
439   for (int i = startView; i <= endView; i++) {
440     is.seekg (sizeHeader + i * sizeView);
441     is.read (reinterpret_cast<char*>(pViewData), sizeView);
442     os.write (reinterpret_cast<char*>(pViewData), sizeView);
443     if (is.fail() || os.fail())
444       break;
445   }
446   
447   delete pViewData;
448   if (is.fail()) 
449     sys_error (ERR_FATAL, "Error reading projection file");
450   if (os.fail()) 
451     sys_error (ERR_FATAL, "Error writing projection file");
452   
453   return (! (is.fail() | os.fail()));
454 }
455
456 bool 
457 Projections::copyHeader (const std::string& filename, std::ostream& os)
458 {
459   return copyHeader (filename.c_str(), os);
460 }
461
462 bool
463 Projections::copyHeader (const char* const filename, std::ostream& os)
464 {
465   frnetorderstream is (filename, std::ios::in | std::ios::binary);
466   kuint16 sizeHeader, signature;
467   is.readInt16 (sizeHeader);
468   is.readInt16 (signature);
469   is.seekg (0);
470   if (signature != m_signature) {
471     sys_error (ERR_FATAL, "Illegal signature in projection file %s", filename);
472     return false;
473   }
474   
475   unsigned char* pHdrData = new unsigned char [sizeHeader];
476   is.read (reinterpret_cast<char*>(pHdrData), sizeHeader);
477   if (is.fail()) {
478     sys_error (ERR_FATAL, "Error reading header");
479     return false;
480   }
481   
482   os.write (reinterpret_cast<char*>(pHdrData), sizeHeader);
483   if (os.fail()) {
484     sys_error (ERR_FATAL, "Error writing header");
485     return false;
486   }
487   
488   return true;
489 }
490
491 bool
492 Projections::write (const std::string& filename)
493 {
494   return write (filename.c_str());
495 }
496
497 bool
498 Projections::write (const char* filename)
499 {
500   frnetorderstream fs (filename, std::ios::out | std::ios::binary | std::ios::trunc | std::ios::ate);
501   m_filename = filename;
502   if (! fs) {
503     sys_error (ERR_SEVERE, "Error opening file %s for output [projections_create]", filename);
504     return false;
505   }
506   
507   if (! headerWrite (fs))
508     return false;
509   
510   if (m_projData != NULL) {
511     for (int i = 0; i < m_nView; i++) {
512       if (! detarrayWrite (fs, *m_projData[i], i))
513         break;
514     }
515   }
516   if (! fs)
517     return false;
518   
519   fs.close();
520   
521   return true;
522 }
523
524 /* NAME
525 *   detarrayRead                Read a Detector Array structure from the disk
526 *
527 * SYNOPSIS
528 *   detarrayRead (proj, darray, view_num)
529 *   DETARRAY *darray            Detector array storage location to be filled
530 *   int      view_num           View number to read
531 */
532
533 bool
534 Projections::detarrayRead (fnetorderstream& fs, DetectorArray& darray, const int iview)
535 {
536   const int detval_bytes = darray.nDet() * sizeof(kfloat32);
537   const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
538   const int view_bytes = detheader_bytes + detval_bytes;
539   const off_t start_data = m_headerSize + (iview * view_bytes);
540   DetectorValue* detval_ptr = darray.detValues();  
541   kfloat64 view_angle;
542   kuint32 nDet;
543   
544   fs.seekg (start_data);
545   
546   fs.readFloat64 (view_angle);
547   fs.readInt32 (nDet);
548   darray.setViewAngle (view_angle);
549   //  darray.setNDet ( nDet);
550   
551   for (unsigned int i = 0; i < nDet; i++) {
552     kfloat32 detval;
553     fs.readFloat32 (detval);
554     detval_ptr[i] = detval;
555   }
556   if (! fs)
557     return false;
558   
559   return true;
560 }
561
562
563 /* NAME
564 *   detarrayWrite                       Write detector array data to the disk
565 *
566 * SYNOPSIS
567 *   detarrayWrite (darray, view_num)
568 *   DETARRAY *darray                    Detector array data to be written
569 *   int      view_num                   View number to write
570 *
571 * DESCRIPTION
572 *       This routine writes the detarray data from the disk sequentially to
573 *    the file that was opened with open_projections().  Data is written in
574 *    binary format.
575 */
576
577 bool
578 Projections::detarrayWrite (fnetorderstream& fs, const DetectorArray& darray, const int iview)
579 {
580   const int detval_bytes = darray.nDet() * sizeof(float);
581   const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
582   const int view_bytes = detheader_bytes + detval_bytes;
583   const off_t start_data = m_headerSize + (iview * view_bytes);
584   const DetectorValue* const detval_ptr = darray.detValues();  
585   kfloat64 view_angle = darray.viewAngle();
586   kuint32 nDet = darray.nDet();
587   
588   fs.seekp (start_data);
589   if (! fs) {
590     sys_error (ERR_SEVERE, "Error seeking detectory array [detarrayWrite]");
591     return false;
592   }
593   
594   fs.writeFloat64 (view_angle);
595   fs.writeInt32 (nDet);
596   
597   for (unsigned int i = 0; i < nDet; i++) {
598     kfloat32 detval = detval_ptr[i];
599     fs.writeFloat32 (detval);
600   }
601   
602   if (! fs)
603     return (false);
604   
605   return true;
606 }
607
608 /* NAME
609 *   printProjectionData                 Print projections data
610 *
611 * SYNOPSIS
612 *   printProjectionData ()
613 */
614
615 void
616 Projections::printProjectionData ()
617 {
618   printProjectionData (0, nView() - 1);
619 }
620
621 void
622 Projections::printProjectionData (int startView, int endView)
623 {
624   printf("Projections Data\n\n");
625   printf("Description: %s\n", m_remark.c_str());
626   printf("Geometry: %s\n", Scanner::convertGeometryIDToName (m_geometry));
627   printf("nView       = %8d           nDet = %8d\n", m_nView, m_nDet);
628   printf("focalLength = %8.4f  fieldOfView = %8.4f\n", m_focalLength, m_fieldOfView);
629   printf("rotStart    = %8.4f       rotInc = %8.4f\n", m_rotStart, m_rotInc);
630   printf("detStart    = %8.4f       detInc = %8.4f\n", m_detStart, m_detInc);
631   if (m_projData != NULL) {
632     if (startView < 0)
633       startView = 0;
634     if (startView > m_nView - 1)
635       startView = m_nView - 1;
636     if (endView > m_nView - 1)
637       endView = m_nView - 1;
638     for (int ir = startView; ir <= endView - 1; ir++) {
639       printf("View %d: angle %f\n", ir, m_projData[ir]->viewAngle());
640       DetectorValue* detval = m_projData[ir]->detValues();
641       for (int id = 0; id < m_projData[ir]->nDet(); id++)
642         printf("%8.4f  ", detval[id]);
643       printf("\n");
644     }
645   }
646 }
647
648 void 
649 Projections::printScanInfo (std::ostringstream& os) const
650 {
651   os << "Number of detectors: " << m_nDet << "\n";
652   os << "    Number of views: " << m_nView<< "\n";
653   os << "             Remark: " << m_remark.c_str()<< "\n";
654   os << "           Geometry: " << Scanner::convertGeometryIDToName (m_geometry)<< "\n";
655   os << "       Focal Length: " << m_focalLength<< "\n";
656   os << "      Field Of View: " << m_fieldOfView<< "\n";
657   os << "             phmLen: " << m_phmLen<< "\n";
658   os << "           detStart: " << m_detStart<< "\n";
659   os << "             detInc: " << m_detInc<< "\n";
660   os << "           rotStart: " << m_rotStart<< "\n";
661   os << "             rotInc: " << m_rotInc<< "\n";
662 }
663
664
665 bool 
666 Projections::convertPolar (ImageFile& rIF, int iInterpolationID)
667 {
668   unsigned int nx = rIF.nx();
669   unsigned int ny = rIF.ny();
670   ImageFileArray v = rIF.getArray();
671   ImageFileArray vImag = rIF.getImaginaryArray();
672
673   if (! v || nx == 0 || ny == 0)
674     return false;
675   
676   Array2d<double> adView (nx, ny);
677   Array2d<double> adDet (nx, ny);
678   double** ppdView = adView.getArray();
679   double** ppdDet = adDet.getArray();
680
681   calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet);
682
683   std::complex<double>** ppcDetValue = new std::complex<double>* [m_nView];
684   for (unsigned int iView = 0; iView < m_nView; iView++) {
685     ppcDetValue[iView] = new std::complex<double> [m_nDet];
686     for (unsigned int iDet = 0; iDet < m_nDet; iDet++)
687       ppcDetValue[iView][iDet] = std::complex<double>(getDetectorArray (iView).detValues()[iDet], 0);
688   }
689
690   interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, m_nView, m_nDet, iInterpolationID);
691
692   for (iView = 0; iView < m_nView; iView++)
693     delete [] ppcDetValue[iView];
694   delete [] ppcDetValue;
695
696   return true;
697 }
698
699
700 bool 
701 Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad)
702 {
703   unsigned int nx = rIF.nx();
704   unsigned int ny = rIF.ny();
705   ImageFileArray v = rIF.getArray();
706   if (! rIF.isComplex())
707     rIF.convertRealToComplex();
708   ImageFileArray vImag = rIF.getImaginaryArray();
709
710   if (! v || nx == 0 || ny == 0)
711     return false;
712   
713   Array2d<double> adView (nx, ny);
714   Array2d<double> adDet (nx, ny);
715   double** ppdView = adView.getArray();
716   double** ppdDet = adDet.getArray();
717
718   std::complex<double>** ppcDetValue = new std::complex<double>* [m_nView];
719   unsigned int iView;
720   double* pdDet = new double [m_nDet];
721   fftw_complex* pcIn = new fftw_complex [m_nDet];
722   fftw_plan plan = fftw_create_plan (m_nDet, FFTW_FORWARD, FFTW_IN_PLACE);
723
724   for (iView = 0; iView < m_nView; iView++) {
725     ppcDetValue[iView] = new std::complex<double> [m_nDet];
726     unsigned int iDet;
727     for (iDet = 0; iDet < m_nDet; iDet++) {
728       pcIn[iDet].re = getDetectorArray(iView).detValues()[iDet];
729       pcIn[iDet].im = 0;
730     }
731     fftw_one (plan, pcIn, NULL);
732     for (iDet = 0; iDet < m_nDet; iDet++)
733       ppcDetValue[iView][iDet] = std::complex<double> (pcIn[iDet].re, pcIn[iDet].im);
734    Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], m_nDet);
735   }
736
737   fftw_destroy_plan (plan);  
738   delete [] pcIn;
739   
740   calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet);
741
742   interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, m_nView, m_nDet, iInterpolationID);
743
744   for (iView = 0; iView < m_nView; iView++)
745     delete [] ppcDetValue[iView];
746   delete [] ppcDetValue;
747
748   return true;
749 }
750
751 void
752 Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double** ppdView, double** ppdDet)
753 {
754   double xMin = -m_phmLen / 2;
755   double xMax = xMin + m_phmLen;
756   double yMin = -m_phmLen / 2;
757   double yMax = yMin + m_phmLen;
758   
759   double xInc = (xMax - xMin) / nx;     // size of cells
760   double yInc = (yMax - yMin) / ny;
761   
762   int iDetCenter = (m_nDet - 1) / 2;    // index refering to L=0 projection 
763
764   if (m_geometry != Scanner::GEOMETRY_PARALLEL) {
765     sys_error (ERR_WARNING, "convertPolar supports Parallel only");
766   }
767   
768   double x = xMin + xInc / 2;   // Rectang coords of center of pixel 
769   for (unsigned int ix = 0; ix < nx; x += xInc, ix++) {
770     double y = yMin + yInc / 2;
771     for (unsigned int iy = 0; iy < ny; y += yInc, iy++) {
772       double r = ::sqrt (x * x + y * y);
773       double phi = atan2 (y, x);
774
775       if (phi >= PI) {
776         phi -= PI;
777       } else if (phi < 0) {
778         phi += PI;
779       } else
780         r = -r;
781       
782       ppdView[ix][iy] = (phi - m_rotStart) / m_rotInc;
783       ppdDet[ix][iy] = (r / m_detInc) + iDetCenter;
784     }
785   }
786 }
787
788 void
789 Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag,
790      unsigned int nx, unsigned int ny, std::complex<double>** ppcDetValue,
791      double** ppdView, double** ppdDet, unsigned int nView, unsigned int nDet, int iInterpolationID)
792 {
793   for (unsigned int ix = 0; ix < ny; ix++) {
794     for (unsigned int iy = 0; iy < ny; iy++) {
795       if (iInterpolationID == POLAR_INTERP_NEAREST) {
796         int iView = nearest<int> (ppdView[ix][iy]);
797         int iDet = nearest<int> (ppdDet[ix][iy]);
798         if (iView == nView) {
799           iView = 0;
800        //   iDet = m_nDet - iDet;
801         }
802         if (iDet >= 0 && iDet < nDet && iView >= 0 && iView < nView) {
803           v[ix][iy] = ppcDetValue[iView][iDet].real();
804           if (vImag)
805             vImag[ix][iy] = ppcDetValue[iView][iDet].imag();
806         } else {
807           sys_error (ERR_SEVERE, "Can't find projection data for ix=%d,iy=%d with radView=%f and radDet=%f", 
808             ix, iy, ppdView[ix][iy], ppdDet[ix][iy]);
809           v[ix][iy] = 0;
810         }
811       } else if (iInterpolationID == POLAR_INTERP_BILINEAR) {
812         int iFloorView = ppdView[ix][iy];
813         double dFracView = ppdView[ix][iy] - iFloorView;
814         int iFloorDet = ppdDet[ix][iy];
815         double dFracDet = ppdDet[ix][iy] - iFloorDet;
816
817         if (iFloorDet >= 0 && iFloorView >= 0) { 
818           std::complex<double> v1 = ppcDetValue[iFloorView][iFloorDet];
819           std::complex<double> v2, v3, v4;
820           if (iFloorView < nView - 1)
821             v2 = ppcDetValue[iFloorView + 1][iFloorDet];
822           else 
823             v2 = ppcDetValue[0][iFloorDet];
824           if (iFloorDet < nDet - 1) 
825             v4 = ppcDetValue[iFloorView][iFloorDet+1];
826           else
827             v4 = v1;
828           if (iFloorView < nView - 1 && iFloorDet < nDet - 1)
829             v3 = ppcDetValue [iFloorView+1][iFloorDet+1];
830           else if (iFloorView < nView - 1)
831             v3 = v2;
832           else 
833             v3 = ppcDetValue[0][iFloorDet+1];
834           std::complex<double> vInterp = (1 - dFracView) * (1 - dFracDet) * v1 +
835             dFracView * (1 - dFracDet) * v2 + dFracView * dFracDet * v3 +
836             dFracDet * (1 - dFracView) * v4;
837           v[ix][iy] = vInterp.real();
838           if (vImag)
839             vImag[ix][iy] = vInterp.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           if (vImag)
845             vImag[ix][iy] = 0;
846         }
847       } else if (iInterpolationID == POLAR_INTERP_BICUBIC) {
848         v[ix][iy] =0;
849           if (vImag)
850             vImag[ix][iy] = 0;
851       }
852     }
853   }
854 }
855
856