r5231: Auto commit for Debian build
[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.83 2003/07/04 21:39:40 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 #include <ctime>
29 #include "interpolator.h"
30
31 const kuint16 Projections::m_signature = ('P'*256 + 'J');
32
33 const int Projections::POLAR_INTERP_INVALID = -1;
34 const int Projections::POLAR_INTERP_NEAREST = 0;
35 const int Projections::POLAR_INTERP_BILINEAR = 1;
36 const int Projections::POLAR_INTERP_BICUBIC = 2;
37
38 const char* const Projections::s_aszInterpName[] = 
39 {
40   "nearest",
41   "bilinear",
42 //  {"bicubic"},
43 };
44
45 const char* const Projections::s_aszInterpTitle[] = 
46 {
47   "Nearest",
48   "Bilinear",
49 //  {"Bicubic"},
50 };
51
52 const int Projections::s_iInterpCount = sizeof(s_aszInterpName) / sizeof(char*);
53
54
55
56 /* NAME
57 *    Projections                Constructor for projections matrix storage 
58 *
59 * SYNOPSIS
60 *    proj = projections_create (filename, nView, nDet)
61 *    Projections& proj          Allocated projections structure & matrix
62 *    int nView                  Number of rotated view
63 *    int nDet                   Number of detectors
64 *
65 */
66
67 Projections::Projections (const Scanner& scanner)
68 : m_projData(0)
69 {
70   initFromScanner (scanner);
71 }
72
73
74 Projections::Projections (const int nView, const int nDet)
75 : m_projData(0)
76 {
77   init (nView, nDet);
78 }
79
80 Projections::Projections (void)
81 : m_projData(0)
82 {
83   init (0, 0);
84 }
85
86 Projections::~Projections (void)
87 {
88   deleteProjData();
89 }
90
91 int
92 Projections::convertInterpNameToID (const char* const interpName)
93 {
94   int interpID = POLAR_INTERP_INVALID;
95   
96   for (int i = 0; i < s_iInterpCount; i++)
97     if (strcasecmp (interpName, s_aszInterpName[i]) == 0) {
98       interpID = i;
99       break;
100     }
101     
102     return (interpID);
103 }
104
105 const char*
106 Projections::convertInterpIDToName (const int interpID)
107 {
108   static const char *interpName = "";
109   
110   if (interpID >= 0 && interpID < s_iInterpCount)
111     return (s_aszInterpName[interpID]);
112   
113   return (interpName);
114 }
115
116 const char*
117 Projections::convertInterpIDToTitle (const int interpID)
118 {
119   static const char *interpTitle = "";
120   
121   if (interpID >= 0 && interpID < s_iInterpCount)
122     return (s_aszInterpTitle[interpID]);
123   
124   return (interpTitle);
125 }
126
127
128
129 void
130 Projections::init (const int nView, const int nDet)
131 {
132   m_label.setLabelType (Array2dFileLabel::L_HISTORY);
133   m_nView = nView;
134   m_nDet = nDet;
135   newProjData ();
136   
137   time_t t = time (NULL);
138   tm* lt = localtime (&t);
139   m_year = lt->tm_year;
140   m_month = lt->tm_mon;
141   m_day = lt->tm_mday;
142   m_hour = lt->tm_hour;
143   m_minute = lt->tm_min;
144   m_second = lt->tm_sec;
145 }
146
147 void
148 Projections::initFromScanner (const Scanner& scanner)
149 {
150   m_label.setLabelType (Array2dFileLabel::L_HISTORY);
151   deleteProjData();
152   init (scanner.nView(), scanner.nDet());
153   
154   m_rotInc = scanner.rotInc();
155   m_detInc = scanner.detInc();
156   m_detStart =  scanner.detStart();
157   m_geometry = scanner.geometry();
158   m_dFocalLength = scanner.focalLength();
159   m_dSourceDetectorLength = scanner.sourceDetectorLength();
160   m_dViewDiameter = scanner.viewDiameter();
161   m_rotStart = scanner.offsetView()*scanner.rotInc();
162   m_dFanBeamAngle = scanner.fanBeamAngle();
163 }
164
165 void
166 Projections::setNView (int nView)  // used by MPI to reduce # of views
167 {
168   deleteProjData();
169   init (nView, m_nDet);
170 }
171
172 //  Helical 180 Linear Interpolation.
173 //  This member function takes a set of helical scan projections and 
174 //  performs a linear interpolation between pairs of complementary rays 
175 //  to produce a single projection data set approximating what would be
176 //  measured at a single axial plane.
177 //  Complementary rays are rays which traverse the same path through the 
178 //  phantom in opposite directions.
179 //
180 //  For parallel beam geometry, a ray with a given gantry angle beta and a
181 //  detector iDet will have a complementary ray at beta + pi and nDet-iDet
182 //
183 //  For equiangular or equilinear beam geometry the complementary ray to
184 //  gantry angle beta and fan-beam angle gamma is at 
185 //  beta-hat = beta +2*gamma + pi, and gamma-hat =  -gamma.
186 //  Note that beta-hat - beta depends on gamma and is not constant.
187 //
188 //  The algorithm used here is from Crawford and King, Med. Phys. 17(6)
189 //  1990 p967; what they called method "C", CSH-HH.  It uses interpolation only
190 //  between pairs of complementary rays on either side of an image plane.
191 //  Input data must sample gantry angles from zero to  
192 //  (2*pi + 2* fan-beam-angle).  The data set produced contains gantry
193 //  angles from 0 to Pi+fan-beam-angle.  This is a "halfscan" data set,
194 //  which still contains redundant data, and can be used with a half scan 
195 //  reconstruction to produce an image.
196 //  In this particular implementation a lower triangle from (beta,gamma) =
197 //  (0,-fanAngle/2)->(2*fanAngle,-fanAngle/2)->(0,fanAngle/2) contains
198 //  zeros, but is actually redundant with data contained in the region
199 //  (pi+fanAngle,-fanAngle/2)->(pi+fanAngle, fanAngle/2) ->(pi-fanAngle,
200 //  fanAngle/2).  
201 //
202 int 
203 Projections::Helical180LI(int interpolation_view)
204 {
205    if (m_geometry == Scanner::GEOMETRY_INVALID) 
206    {
207        std::cerr << "Invalid geometry " << m_geometry << std::endl;
208        return (2);
209    } 
210    else if (m_geometry == Scanner::GEOMETRY_PARALLEL) 
211    {
212        std::cerr << "Helical 180LI not yet implemented for PARALLEL geometry"
213                    << std::endl;
214        return (2);
215    }
216    else if (m_geometry == Scanner::GEOMETRY_EQUILINEAR) 
217    {
218        std::cerr << "Helical 180LI not yet implemented for EQUILINEAR geometry"
219                    << std::endl;
220        return (2);
221    }
222    else if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR) 
223    {
224            return Helical180LI_Equiangular(interpolation_view);
225    }
226    else
227    {
228        std::cerr << "Invalid geometry  in projection data file" << m_geometry
229                    << std::endl;
230        return (2);
231    }
232 }
233 int
234 Projections::Helical180LI_Equiangular(int interpView)
235 {
236    double dbeta = m_rotInc; 
237    double dgamma =  m_detInc; 
238    double fanAngle = m_dFanBeamAngle;
239    int offsetView=0;
240    
241    // is there enough data in the data set?  Should have 2(Pi+fanAngle)
242    // coverage minimum
243    if ( m_nView <  static_cast<int>((2*( PI + fanAngle ) ) / dbeta) -1 ){
244        std::cerr   << "Data set does not include 360 +2*FanBeamAngle views"
245                    << std::endl;
246        return (1);
247    }
248
249    if (interpView < 0)   // use default position at PI+fanAngle
250    {
251        interpView = static_cast<int> ((PI+fanAngle)/dbeta);
252    }
253    else
254    {
255        // check if there is PI+fanAngle data on either side of the 
256        // of the specified image plane
257        if ( interpView*dbeta < PI+fanAngle ||            
258             interpView*dbeta + PI + fanAngle > m_nView*dbeta) 
259        {
260            std::cerr << "There isn't PI+fanAngle of data on either side of the requested interpolation view" << std::endl;
261            return(1);
262        }
263        offsetView = interpView - static_cast<int>((PI+fanAngle)/dbeta);
264
265    }
266    int last_interp_view = static_cast<int> ((PI+fanAngle)/dbeta);
267
268    
269 // make a new array for data...
270    class DetectorArray ** newdetarray = new DetectorArray * [last_interp_view+1];
271    for ( int i=0 ; i <= last_interp_view ; i++ ){
272        newdetarray[i] = new DetectorArray (m_nDet);
273        newdetarray[i]->setViewAngle((i+offsetView)*dbeta);
274        DetectorValue* newdetval = (newdetarray[i])->detValues();
275        // and initialize the data to zero
276        for (int j=0; j < m_nDet; j++) 
277            newdetval[j] = 0.;
278    }
279
280    int last_acq_view = 2*last_interp_view;
281    for ( int iView = 0 ; iView <= last_acq_view; iView++) {
282        double beta = iView * dbeta; 
283        
284        for ( int iDet = 0; iDet < m_nDet; iDet++) {
285            double gamma = (iDet -(m_nDet-1)/2)* dgamma ;
286            int newiView, newiDet;
287            if (beta < PI+fanAngle) { //if (PI +fanAngle - beta > dbeta )  
288                //newbeta = beta; 
289                //newgamma = gamma; 
290                newiDet = iDet; 
291                newiView = iView; 
292            }
293            else // (beta > PI+fanAngle)
294            {
295                //newbeta = beta +2*gamma - 180;
296                //newgamma = -gamma;
297                newiDet = -iDet + (m_nDet -1);
298                // newiView = nearest<int>((beta + 2*gamma - PI)/dbeta);
299                //newiView = static_cast<int>(( (iView*dbeta) + 2*(iDet-(m_nDet-1)/2)*dgamma - PI)/dbeta);
300                newiView = nearest<int>(( (iView*dbeta) + 2*(iDet-(m_nDet-1)/2)*dgamma - PI)/dbeta);
301            } 
302
303 #ifdef DEBUG
304 //std::cout << beta << " "<< gamma << " " << newbeta << " " << newgamma <<"    " << iView-offsetView << " " << iDet << " " << newiView << " " << newiDet << std::endl;
305 //std::cout << iView-offsetView << " " << iDet << " " << newiView << " " << newiDet << std::endl;
306 #endif
307
308            if (   ( beta > fanAngle - 2*gamma) 
309                && ( beta < 2*PI + fanAngle -2*gamma)  )
310           {  // not in region  1 or 8
311                DetectorValue* detval = (m_projData[iView+offsetView])->detValues();
312                DetectorValue* newdetval = (newdetarray[newiView])->detValues();
313                if (   beta > fanAngle - 2*gamma  
314                    && beta <= 2*fanAngle ) {  // in region 2
315                    newdetval[newiDet] += 
316                        (beta +2*gamma - fanAngle)/(PI+2*gamma)
317                                * detval[iDet];
318                } else if ( beta > 2*fanAngle  
319                           && beta <= PI - 2*gamma) {  // in region 3
320                    newdetval[newiDet] += 
321                        (beta +2*gamma - fanAngle)/(PI+2*gamma)
322                                * detval[iDet];
323                } 
324                else if (   beta > PI -2*gamma  
325                         && beta <= PI + fanAngle ) {  // in region 4
326                    newdetval[newiDet] += 
327                        (beta +2*gamma - fanAngle)/(PI+2*gamma)
328                                * detval[iDet];
329                } 
330                else if (   beta > PI + fanAngle  
331                         && beta <= PI +2*fanAngle -2*gamma) { // in region 5
332                    newdetval[newiDet] += 
333                        (2*PI - beta - 2*gamma + fanAngle)/(PI-2*gamma)
334                                *detval[iDet];
335                } 
336                else if (   beta > PI +2*fanAngle -2*gamma 
337                         && beta <= 2*PI) {  // in region 6
338                    newdetval[newiDet] += 
339                        (2*PI - beta - 2*gamma + fanAngle)/(PI-2*gamma)
340                        *detval[iDet];
341                } 
342                else if (   beta > 2*PI 
343                         && beta <= 2*PI + fanAngle -2*gamma){ // in region 7
344                    newdetval[newiDet] += 
345                        (2*PI - beta -2*gamma + fanAngle)/(PI-2*gamma)
346                        *detval[iDet];
347                } 
348                else 
349                { 
350                    ; // outside region of interest
351                }
352            }
353        }
354    }
355    deleteProjData();
356    m_projData = newdetarray;
357    m_nView = last_interp_view+1;
358
359    return (0); 
360 }
361 // HalfScanFeather:
362 // A HalfScan Projection Data Set  for equiangular geometry, 
363 // covering gantry angles from 0 to  pi+fanBeamAngle 
364 // and fan angle gamma from -fanBeamAngle/2 to fanBeamAngle/2
365 // contains redundant information.  If one copy of this data is left as
366 // zero, (as in the Helical180LI routine above) overweighting is avoided, 
367 // but the discontinuity in the data introduces ringing in the image. 
368 // This routine makes a copy of the data and applies a weighting to avoid
369 // over-representation, as given in Appendix C of Crawford and King, Med
370 // Phys 17 1990, p967.
371 int
372 Projections::HalfScanFeather(void)
373 {
374    double dbeta = m_rotInc; 
375    double dgamma =  m_detInc; 
376    double fanAngle = m_dFanBeamAngle;
377
378 // is there enough data?  
379    if ( m_nView !=  static_cast<int>(( PI+fanAngle ) / dbeta) +1 ){
380        std::cerr   << "Data set does seem have enough data to be a halfscan data set"  << std::endl;
381        return (1);
382    }
383    if (m_geometry == Scanner::GEOMETRY_INVALID) {
384        std::cerr << "Invalid geometry " << m_geometry << std::endl;
385        return (2);
386    }
387
388    if (m_geometry == Scanner::GEOMETRY_PARALLEL) {
389        std::cerr << "HalfScanFeather not yet implemented for PARALLEL geometry"<< std::endl;
390        return (2);
391    }
392
393    for ( int iView2 = 0 ; iView2 < m_nView; iView2++) {
394        double beta2 = iView2 * dbeta; 
395        for ( int iDet2 = 0; iDet2 < m_nDet; iDet2++) {
396            double gamma2 = (iDet2 -(m_nDet-1)/2)* dgamma ;
397            if ( ( beta2 >= PI  - 2*gamma2) ) {  // in redundant data region 
398                int iView1, iDet1;
399                iDet1 =  (m_nDet -1) - iDet2;
400                //iView1 = nearest<int>((beta2 + 2*gamma2 - PI)/dbeta);
401                iView1 = nearest<int>(( (iView2*dbeta) 
402                                + 2*(iDet2-(m_nDet-1)/2)*dgamma - PI)/dbeta);
403
404
405                DetectorValue* detval2 = (m_projData[iView2])->detValues();
406                DetectorValue* detval1 = (m_projData[iView1])->detValues();
407
408                detval1[iDet1] = detval2[iDet2] ;
409
410                double x, w1,w2,beta1, gamma1;
411                beta1= iView1*dbeta; 
412                gamma1 = -gamma2;
413                if ( beta1 <= (fanAngle - 2*gamma1) )
414                    x = beta1 / ( fanAngle - 2*gamma1);
415                else if ( (fanAngle  - 2*gamma1 <= beta1 ) && beta1 <= PI - 2*gamma1) 
416                    x = 1; 
417                else if ( (PI - 2*gamma1 <= beta1 ) && ( beta1 <=PI + fanAngle) )  
418                    x = (PI +fanAngle - beta1)/(fanAngle + 2*gamma1);
419                else {
420                    std::cerr << "Shouldn't be here!"<< std::endl;
421                    return(4);
422                }
423                w1 = (3*x - 2*x*x)*x;
424                w2 = 1-w1;
425                detval1[iDet1] *= w1; 
426                detval2[iDet2] *= w2;
427
428            } 
429        }
430    }
431    // heuristic scaling, why this factor?  
432    double scalefactor = m_nView * m_rotInc / PI;
433    for ( int iView = 0 ; iView < m_nView; iView++) {
434        DetectorValue* detval = (m_projData[iView])->detValues();
435        for ( int iDet = 0; iDet < m_nDet; iDet++) {
436            detval[iDet] *= scalefactor;
437        }
438    }
439
440    return (0); 
441 }
442
443 // NAME
444 // newProjData
445
446 void 
447 Projections::newProjData (void)
448 {
449   if (m_projData)
450     sys_error(ERR_WARNING, "m_projData != NULL [newProjData]");
451   
452   if (m_nView > 0 && m_nDet) {
453     m_projData = new DetectorArray* [m_nView];
454     
455     for (int i = 0; i < m_nView; i++)
456       m_projData[i] = new DetectorArray (m_nDet);
457   }
458 }
459
460
461 /* NAME
462 *   projections_free                    Free memory allocated to projections
463 *
464 * SYNOPSIS
465 *   projections_free(proj)
466 *   Projections& proj                           Projectionss to be deallocated
467 */
468
469 void 
470 Projections::deleteProjData (void)
471 {
472   if (m_projData != NULL) {
473     for (int i = 0; i < m_nView; i++)
474       delete m_projData[i];
475     
476     delete m_projData;
477     m_projData = NULL;
478   }
479 }
480
481
482 /* NAME
483 *    Projections::headerWwrite         Write data header for projections file
484 *
485 */
486
487 bool 
488 Projections::headerWrite (fnetorderstream& fs)
489 {
490   kuint16 _hsize = m_headerSize;
491   kuint16 _signature = m_signature;
492   kuint32 _nView = m_nView;
493   kuint32 _nDet = m_nDet;
494   kuint32 _geom = m_geometry;
495   kuint16 _remarksize = m_remark.length();
496   kuint16 _year = m_year;
497   kuint16 _month = m_month;
498   kuint16 _day = m_day;
499   kuint16 _hour = m_hour;
500   kuint16 _minute = m_minute;
501   kuint16 _second = m_second;
502   
503   kfloat64 _calcTime = m_calcTime;
504   kfloat64 _rotStart = m_rotStart;
505   kfloat64 _rotInc = m_rotInc;
506   kfloat64 _detStart = m_detStart;
507   kfloat64 _detInc = m_detInc;
508   kfloat64 _viewDiameter = m_dViewDiameter;
509   kfloat64 _focalLength = m_dFocalLength;
510   kfloat64 _sourceDetectorLength = m_dSourceDetectorLength;
511   kfloat64 _fanBeamAngle = m_dFanBeamAngle;
512
513   fs.seekp(0);
514   if (! fs)
515     return false;
516   
517   fs.writeInt16 (_hsize);
518   fs.writeInt16 (_signature);
519   fs.writeInt32 (_nView);
520   fs.writeInt32 (_nDet);
521   fs.writeInt32 (_geom);
522   fs.writeFloat64 (_calcTime);
523   fs.writeFloat64 (_rotStart);
524   fs.writeFloat64 (_rotInc);
525   fs.writeFloat64 (_detStart);
526   fs.writeFloat64 (_detInc);
527   fs.writeFloat64 (_viewDiameter);
528   fs.writeFloat64 (_focalLength);
529   fs.writeFloat64 (_sourceDetectorLength);
530   fs.writeFloat64 (_fanBeamAngle);
531   fs.writeInt16 (_year);
532   fs.writeInt16 (_month);
533   fs.writeInt16 (_day);
534   fs.writeInt16 (_hour);
535   fs.writeInt16 (_minute);
536   fs.writeInt16 (_second);
537   fs.writeInt16 (_remarksize);
538   fs.write (m_remark.c_str(), _remarksize);
539   
540   m_headerSize = fs.tellp();
541   _hsize = m_headerSize;
542   fs.seekp(0);
543   fs.writeInt16 (_hsize);
544   if (! fs)
545     return false;
546   
547   return true;
548 }
549
550 /* NAME
551 *    projections_read_header         Read data header for projections file
552 *
553 */
554 bool
555 Projections::headerRead (fnetorderstream& fs)
556 {
557   kuint16 _hsize, _signature, _year, _month, _day, _hour, _minute, _second, _remarksize = 0;
558   kuint32 _nView, _nDet, _geom;
559   kfloat64 _calcTime, _rotStart, _rotInc, _detStart, _detInc, _focalLength, _sourceDetectorLength, _viewDiameter, _fanBeamAngle;
560   
561   fs.seekg(0);
562   if (! fs)
563     return false;
564   
565   fs.readInt16 (_hsize);
566   fs.readInt16 (_signature);
567   fs.readInt32 (_nView);
568   fs.readInt32 (_nDet);
569   fs.readInt32 (_geom);
570   fs.readFloat64 (_calcTime);
571   fs.readFloat64 (_rotStart);
572   fs.readFloat64 (_rotInc);
573   fs.readFloat64 (_detStart);
574   fs.readFloat64 (_detInc);
575   fs.readFloat64 (_viewDiameter);
576   fs.readFloat64 (_focalLength);
577   fs.readFloat64 (_sourceDetectorLength);
578   fs.readFloat64 (_fanBeamAngle);
579   fs.readInt16 (_year);
580   fs.readInt16 (_month);
581   fs.readInt16 (_day);
582   fs.readInt16 (_hour);
583   fs.readInt16 (_minute);
584   fs.readInt16 (_second);
585   fs.readInt16 (_remarksize);
586   
587   if (! fs) {
588     sys_error (ERR_SEVERE, "Error reading header information , _remarksize=%d [projections_read_header]", _remarksize);
589     return false;
590   }
591   
592   if (_signature != m_signature) {
593     sys_error (ERR_SEVERE, "File %s does not have a valid projection file signature", m_filename.c_str());
594     return false;
595   }
596   
597   char* pszRemarkStorage = new char [_remarksize+1];
598   fs.read (pszRemarkStorage, _remarksize);
599   if (! fs) {
600     sys_error (ERR_SEVERE, "Error reading remark, _remarksize = %d", _remarksize);
601     return false;
602   }
603   pszRemarkStorage[_remarksize] = 0;
604   m_remark = pszRemarkStorage;
605   delete pszRemarkStorage;
606   
607   off_t _hsizeread = fs.tellg();
608   if (!fs || _hsizeread != _hsize) {
609     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);
610     return false;
611   }
612   
613   m_headerSize = _hsize;
614   m_nView = _nView;
615   m_nDet = _nDet;
616   m_geometry = _geom;
617   m_calcTime = _calcTime;
618   m_rotStart = _rotStart;
619   m_rotInc = _rotInc;
620   m_detStart = _detStart;
621   m_detInc = _detInc;
622   m_dFocalLength = _focalLength;
623   m_dSourceDetectorLength = _sourceDetectorLength;
624   m_dViewDiameter = _viewDiameter;
625   m_dFanBeamAngle = _fanBeamAngle;
626   m_year = _year;
627   m_month = _month;
628   m_day = _day;
629   m_hour = _hour;
630   m_minute = _minute;
631   m_second = _second;
632   
633   m_label.setLabelType (Array2dFileLabel::L_HISTORY);
634   m_label.setLabelString (m_remark);
635   m_label.setCalcTime (m_calcTime);
636   m_label.setDateTime (m_year, m_month, m_day, m_hour, m_minute, m_second);
637   
638   return true;
639 }
640
641 bool
642 Projections::read (const std::string& filename)
643 {
644   return read (filename.c_str());
645 }
646
647
648 bool
649 Projections::read (const char* filename)
650 {
651   m_filename = filename;
652 #ifdef MSVC
653   frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary);
654 #else
655   frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary); // | std::ios::nocreate);
656 #endif
657   
658   if (fileRead.fail())
659     return false;
660   
661   if (! headerRead (fileRead))
662     return false;
663   
664   deleteProjData ();
665   newProjData();
666   
667   for (int i = 0; i < m_nView; i++) {
668     if (! detarrayRead (fileRead, *m_projData[i], i))
669       break;
670   }
671   
672   fileRead.close();
673   return true;
674 }
675
676
677 bool 
678 Projections::copyViewData (const std::string& filename, std::ostream& os, int startView, int endView)
679 {
680   return copyViewData (filename.c_str(), os, startView, endView);
681 }
682
683 bool 
684 Projections::copyViewData (const char* const filename, std::ostream& os, int startView, int endView)
685 {
686   frnetorderstream is (filename, std::ios::in | std::ios::binary);
687   kuint16 sizeHeader, signature;
688   kuint32 _nView, _nDet;
689   
690   is.seekg (0);
691   if (is.fail()) {
692     sys_error (ERR_SEVERE, "Unable to read projection file %s", filename);
693     return false;
694   }
695
696   is.readInt16 (sizeHeader);
697   is.readInt16 (signature);
698   is.readInt32 (_nView);
699   is.readInt32 (_nDet);
700   int nView = _nView;
701   int nDet = _nDet;
702   
703   if (signature != m_signature) {
704     sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
705     return false;
706   }
707   
708   if (startView < 0)
709     startView = 0;
710   if (startView > nView - 1)
711     startView = nView;
712   if (endView < 0 || endView > nView - 1)
713     endView = nView - 1;
714   
715   if (startView > endView) { // swap if start > end
716     int tempView = endView;
717     endView = startView;
718     startView = tempView;
719   }
720   
721   int sizeView = 8 /* view_angle */ + 4 /* nDet */ + (4 * nDet);
722   unsigned char* pViewData = new unsigned char [sizeView];
723   
724   for (int i = startView; i <= endView; i++) {
725     is.seekg (sizeHeader + i * sizeView);
726     is.read (reinterpret_cast<char*>(pViewData), sizeView);
727     os.write (reinterpret_cast<char*>(pViewData), sizeView);
728     if (is.fail() || os.fail())
729       break;
730   }
731   
732   delete pViewData;
733   if (is.fail()) 
734     sys_error (ERR_SEVERE, "Error reading projection file");
735   if (os.fail()) 
736     sys_error (ERR_SEVERE, "Error writing projection file");
737   
738   return (! (is.fail() | os.fail()));
739 }
740
741 bool 
742 Projections::copyHeader (const std::string& filename, std::ostream& os)
743 {
744   return copyHeader (filename.c_str(), os);
745 }
746
747 bool
748 Projections::copyHeader (const char* const filename, std::ostream& os)
749 {
750   frnetorderstream is (filename, std::ios::in | std::ios::binary);
751   kuint16 sizeHeader, signature;
752   is.readInt16 (sizeHeader);
753   is.readInt16 (signature);
754   is.seekg (0);
755   if (signature != m_signature) {
756     sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
757     return false;
758   }
759   
760   unsigned char* pHdrData = new unsigned char [sizeHeader];
761   is.read (reinterpret_cast<char*>(pHdrData), sizeHeader);
762   if (is.fail()) {
763     sys_error (ERR_SEVERE, "Error reading header");
764     return false;
765   }
766   
767   os.write (reinterpret_cast<char*>(pHdrData), sizeHeader);
768   if (os.fail()) {
769     sys_error (ERR_SEVERE, "Error writing header");
770     return false;
771   }
772   
773   return true;
774 }
775
776 bool
777 Projections::write (const std::string& filename)
778 {
779   return write (filename.c_str());
780 }
781
782 bool
783 Projections::write (const char* filename)
784 {
785   frnetorderstream fs (filename, std::ios::out | std::ios::binary | std::ios::trunc | std::ios::ate);
786   m_filename = filename;
787   if (! fs) {
788     sys_error (ERR_SEVERE, "Error opening file %s for output [projections_create]", filename);
789     return false;
790   }
791   
792   if (! headerWrite (fs))
793     return false;
794   
795   if (m_projData != NULL) {
796     for (int i = 0; i < m_nView; i++) {
797       if (! detarrayWrite (fs, *m_projData[i], i))
798         break;
799     }
800   }
801   if (! fs)
802     return false;
803   
804   fs.close();
805   
806   return true;
807 }
808
809 /* NAME
810 *   detarrayRead                Read a Detector Array structure from the disk
811 *
812 * SYNOPSIS
813 *   detarrayRead (proj, darray, view_num)
814 *   DETARRAY *darray            Detector array storage location to be filled
815 *   int      view_num           View number to read
816 */
817
818 bool
819 Projections::detarrayRead (fnetorderstream& fs, DetectorArray& darray, const int iview)
820 {
821   const int detval_bytes = darray.nDet() * sizeof(kfloat32);
822   const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
823   const int view_bytes = detheader_bytes + detval_bytes;
824   const off_t start_data = m_headerSize + (iview * view_bytes);
825   DetectorValue* detval_ptr = darray.detValues();  
826   kfloat64 view_angle;
827   kuint32 nDet;
828   
829   fs.seekg (start_data);
830   
831   fs.readFloat64 (view_angle);
832   fs.readInt32 (nDet);
833   darray.setViewAngle (view_angle);
834   //  darray.setNDet ( nDet);
835   
836   for (unsigned int i = 0; i < nDet; i++) {
837     kfloat32 detval;
838     fs.readFloat32 (detval);
839     detval_ptr[i] = detval;
840   }
841   if (! fs)
842     return false;
843   
844   return true;
845 }
846
847
848 /* NAME
849 *   detarrayWrite                       Write detector array data to the disk
850 *
851 * SYNOPSIS
852 *   detarrayWrite (darray, view_num)
853 *   DETARRAY *darray                    Detector array data to be written
854 *   int      view_num                   View number to write
855 *
856 * DESCRIPTION
857 *       This routine writes the detarray data from the disk sequentially to
858 *    the file that was opened with open_projections().  Data is written in
859 *    binary format.
860 */
861
862 bool
863 Projections::detarrayWrite (fnetorderstream& fs, const DetectorArray& darray, const int iview)
864 {
865   const int detval_bytes = darray.nDet() * sizeof(float);
866   const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
867   const int view_bytes = detheader_bytes + detval_bytes;
868   const off_t start_data = m_headerSize + (iview * view_bytes);
869   const DetectorValue* const detval_ptr = darray.detValues();  
870   kfloat64 view_angle = darray.viewAngle();
871   kuint32 nDet = darray.nDet();
872   
873   fs.seekp (start_data);
874   if (! fs) {
875     sys_error (ERR_SEVERE, "Error seeking detectory array [detarrayWrite]");
876     return false;
877   }
878   
879   fs.writeFloat64 (view_angle);
880   fs.writeInt32 (nDet);
881   
882   for (unsigned int i = 0; i < nDet; i++) {
883     kfloat32 detval = detval_ptr[i];
884     fs.writeFloat32 (detval);
885   }
886   
887   if (! fs)
888     return (false);
889   
890   return true;
891 }
892
893 /* NAME
894 *   printProjectionData                 Print projections data
895 *
896 * SYNOPSIS
897 *   printProjectionData ()
898 */
899
900 void
901 Projections::printProjectionData ()
902 {
903   printProjectionData (0, nView() - 1);
904 }
905
906 void
907 Projections::printProjectionData (int startView, int endView)
908 {
909   printf("Projections Data\n\n");
910   printf("Description: %s\n", m_remark.c_str());
911   printf("Geometry: %s\n", Scanner::convertGeometryIDToName (m_geometry));
912   printf("nView       = %8d             nDet = %8d\n", m_nView, m_nDet);
913   printf("focalLength = %8.4f   ViewDiameter = %8.4f\n", m_dFocalLength, m_dViewDiameter);
914   printf("fanBeamAngle= %8.4f SourceDetector = %8.4f\n", convertRadiansToDegrees(m_dFanBeamAngle), m_dSourceDetectorLength);
915   printf("rotStart    = %8.4f         rotInc = %8.4f\n", m_rotStart, m_rotInc);
916   printf("detStart    = %8.4f         detInc = %8.4f\n", m_detStart, m_detInc);
917   if (m_projData != NULL) {
918     if (startView < 0)
919       startView = 0;
920     if (endView < 0)
921       endView = m_nView - 1;
922     if (startView > m_nView - 1)
923       startView = m_nView - 1;
924     if (endView > m_nView - 1)
925       endView = m_nView - 1;
926     for (int ir = startView; ir <= endView - 1; ir++) {
927       printf("View %d: angle %f\n", ir, m_projData[ir]->viewAngle());
928       DetectorValue* detval = m_projData[ir]->detValues();
929       for (int id = 0; id < m_projData[ir]->nDet(); id++)
930         printf("%8.4f  ", detval[id]);
931       printf("\n");
932     }
933   }
934 }
935
936 void 
937 Projections::printScanInfo (std::ostringstream& os) const
938 {
939   os << "Number of detectors: " << m_nDet << "\n";
940   os << "Number of views: " << m_nView<< "\n";
941   os << "Description: " << m_remark.c_str()<< "\n";
942   os << "Geometry: " << Scanner::convertGeometryIDToName (m_geometry)<< "\n";
943   os << "Focal Length: " << m_dFocalLength<< "\n";
944   os << "Source Detector Length: " << m_dSourceDetectorLength << "\n";
945   os << "View Diameter: " << m_dViewDiameter<< "\n";
946   os << "Fan Beam Angle: " << convertRadiansToDegrees(m_dFanBeamAngle) << "\n";
947   os << "detStart: " << m_detStart<< "\n";
948   os << "detInc: " << m_detInc<< "\n";
949   os << "rotStart: " << m_rotStart<< "\n";
950   os << "rotInc: " << m_rotInc<< "\n";
951 }
952
953
954 bool 
955 Projections::convertPolar (ImageFile& rIF, int iInterpolationID)
956 {
957   unsigned int nx = rIF.nx();
958   unsigned int ny = rIF.ny();
959   ImageFileArray v = rIF.getArray();
960   ImageFileArray vImag = rIF.getImaginaryArray();
961
962   if (! v || nx == 0 || ny == 0)
963     return false;
964
965   Projections* pProj = this;
966   if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
967     pProj = interpolateToParallel();
968   
969   Array2d<double> adView (nx, ny);
970   Array2d<double> adDet (nx, ny);
971   double** ppdView = adView.getArray();
972   double** ppdDet = adDet.getArray();
973
974   std::complex<double>** ppcDetValue = new std::complex<double>* [pProj->m_nView];
975   int iView;
976   for (iView = 0; iView < pProj->m_nView; iView++) {
977     ppcDetValue[iView] = new std::complex<double> [pProj->m_nDet];
978     DetectorValue* detval = pProj->getDetectorArray (iView).detValues();
979     for (int iDet = 0; iDet < pProj->m_nDet; iDet++)
980       ppcDetValue[iView][iDet] = std::complex<double>(detval[iDet], 0);
981   }
982
983   pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet, pProj->m_nDet, 1., pProj->m_detInc);
984
985   pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet, 
986     pProj->m_nDet, iInterpolationID);
987
988   for (iView = 0; iView < pProj->m_nView; iView++)
989     delete [] ppcDetValue[iView];
990   delete [] ppcDetValue;
991
992   if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
993     delete pProj;
994
995   return true;
996 }
997
998
999 bool 
1000 Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad)
1001 {
1002 #ifndef HAVE_FFTW
1003   rIF.arrayDataClear();
1004   return false;
1005 #else
1006   unsigned int nx = rIF.nx();
1007   unsigned int ny = rIF.ny();
1008   ImageFileArray v = rIF.getArray();
1009   if (! rIF.isComplex())
1010     rIF.convertRealToComplex();
1011   ImageFileArray vImag = rIF.getImaginaryArray();
1012
1013   if (! v || nx == 0 || ny == 0)
1014     return false;
1015   
1016   Projections* pProj = this;
1017   if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
1018     pProj = interpolateToParallel();
1019
1020   int iInterpDet = static_cast<int>(static_cast<double>(sqrt(nx*nx+ny*ny)));
1021   int iNumInterpDetWithZeros = ProcessSignal::addZeropadFactor (iInterpDet, iZeropad);
1022   double dProjScale = iInterpDet / (pProj->viewDiameter() * 0.5);
1023   double dZeropadRatio = static_cast<double>(iNumInterpDetWithZeros) / static_cast<double>(iInterpDet);
1024
1025   fftw_plan plan = fftw_create_plan (iNumInterpDetWithZeros, FFTW_FORWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
1026
1027   fftw_complex* pcIn = new fftw_complex [iNumInterpDetWithZeros];
1028   std::complex<double>** ppcDetValue = new std::complex<double>* [pProj->m_nView];
1029   //double dInterpScale = (pProj->m_nDet-1) / static_cast<double>(iInterpDet-1);
1030   double dInterpScale = pProj->m_nDet / static_cast<double>(iInterpDet);
1031   
1032   double dFFTScale = 1. / static_cast<double>(iInterpDet * iInterpDet);
1033   int iMidPoint = iInterpDet / 2;
1034   double dMidPoint = static_cast<double>(iInterpDet) / 2.;
1035   int iZerosAdded = iNumInterpDetWithZeros - iInterpDet;
1036
1037   // For each view, interpolate, shift to center at origin, and FFT
1038   for (int iView = 0; iView < m_nView; iView++) {
1039     DetectorValue* detval = pProj->getDetectorArray(iView).detValues();
1040     LinearInterpolator<DetectorValue> projInterp (detval, pProj->m_nDet);
1041     for (int iDet = 0; iDet < iInterpDet; iDet++) {
1042       double dInterpPos = (m_nDet / 2.) + (iDet - dMidPoint) * dInterpScale;
1043       pcIn[iDet].re = projInterp.interpolate (dInterpPos) * dProjScale;
1044       pcIn[iDet].im = 0;
1045     }
1046
1047     Fourier::shuffleFourierToNaturalOrder (pcIn, iInterpDet);
1048     if (iZerosAdded > 0) {
1049        for (int iDet1 = iInterpDet -1; iDet1 >= iMidPoint; iDet1--)
1050         pcIn[iDet1+iZerosAdded] = pcIn[iDet1];
1051       for (int iDet2 = iMidPoint; iDet2 < iMidPoint + iZerosAdded; iDet2++) 
1052         pcIn[iDet2].re = pcIn[iDet2].im = 0;
1053     }
1054
1055     fftw_one (plan, pcIn, NULL);
1056
1057     ppcDetValue[iView] = new std::complex<double> [iNumInterpDetWithZeros];
1058     for (int iD = 0; iD < iNumInterpDetWithZeros; iD++) {
1059       ppcDetValue[iView][iD] = std::complex<double> (pcIn[iD].re * dFFTScale, pcIn[iD].im * dFFTScale); 
1060     }
1061
1062     Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], iNumInterpDetWithZeros);
1063   }
1064   delete [] pcIn;
1065
1066   fftw_destroy_plan (plan);  
1067   
1068   Array2d<double> adView (nx, ny);
1069   Array2d<double> adDet (nx, ny);
1070   double** ppdView = adView.getArray();
1071   double** ppdDet = adDet.getArray();
1072   pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet, iNumInterpDetWithZeros, dZeropadRatio, 
1073     pProj->m_detInc * dInterpScale);
1074
1075   pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet, 
1076     iNumInterpDetWithZeros, iInterpolationID);
1077
1078   if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
1079     delete pProj;
1080
1081   for (int i = 0; i < m_nView; i++)
1082     delete [] ppcDetValue[i];
1083   delete [] ppcDetValue;
1084
1085   return true;
1086 #endif
1087 }
1088
1089
1090 void
1091 Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double** ppdView, double** ppdDet,
1092                                         int iNumDetWithZeros, double dZeropadRatio, double dDetInc)
1093 {
1094   double dLength = viewDiameter();
1095   double xMin = -dLength / 2;
1096   double xMax = xMin + dLength;
1097   double yMin = -dLength / 2;
1098   double yMax = yMin + dLength;
1099   double xCent = (xMin + xMax) / 2;
1100   double yCent = (yMin + yMax) / 2;
1101
1102   xMin = (xMin - xCent) * dZeropadRatio + xCent;
1103   xMax = (xMax - xCent) * dZeropadRatio + xCent;
1104   yMin = (yMin - yCent) * dZeropadRatio + yCent;
1105   yMax = (yMax - yCent) * dZeropadRatio + yCent;
1106
1107   double xInc = (xMax - xMin) / nx;     // size of cells
1108   double yInc = (yMax - yMin) / ny;
1109
1110   double dDetCenter = (iNumDetWithZeros - 1) / 2.;      // index refering to L=0 projection 
1111   // +1 is correct for frequency data, ndet-1 is correct for projections
1112   //  if (isEven (iNumDetWithZeros))
1113   //    dDetCenter = (iNumDetWithZeros + 0) / 2;        
1114
1115   // Calculates polar coordinates (view#, det#) for each point on phantom grid
1116   double x = xMin + xInc / 2;   // Rectang coords of center of pixel 
1117   for (unsigned int ix = 0; ix < nx; x += xInc, ix++) {
1118     double y = yMin + yInc / 2;
1119     for (unsigned int iy = 0; iy < ny; y += yInc, iy++) {
1120       double r = ::sqrt (x * x + y * y);
1121       double phi = atan2 (y, x);
1122
1123       if (phi <= -m_rotInc / 2)
1124         phi += TWOPI;
1125       if (phi >= PI - (m_rotInc / 2)) {
1126         phi -= PI;
1127         r = -r;
1128       }
1129       
1130       ppdView[ix][iy] = (phi - m_rotStart) / m_rotInc;
1131       ppdDet[ix][iy] = (r / dDetInc) + dDetCenter;
1132     }
1133   }
1134 }
1135
1136 void
1137 Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag,
1138      unsigned int nx, unsigned int ny, std::complex<double>** ppcDetValue, double** ppdView, 
1139      double** ppdDet, unsigned int nView, unsigned int nDet, unsigned int nDetWithZeros, int iInterpolationID)
1140 {
1141   typedef std::complex<double> complexValue;
1142
1143   BilinearPolarInterpolator<complexValue>* pBilinear = NULL;  
1144   BicubicPolyInterpolator<complexValue>* pBicubic = NULL;  
1145   if (iInterpolationID == POLAR_INTERP_BILINEAR)
1146     pBilinear = new BilinearPolarInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
1147   else if (iInterpolationID == POLAR_INTERP_BICUBIC)
1148     pBicubic = new BicubicPolyInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
1149   
1150   for (unsigned int ix = 0; ix < ny; ix++) {
1151     for (unsigned int iy = 0; iy < ny; iy++) {
1152       if (iInterpolationID == POLAR_INTERP_NEAREST) {
1153         unsigned int iView = nearest<int> (ppdView[ix][iy]);
1154         unsigned int iDet = nearest<int> (ppdDet[ix][iy]);
1155         if (iView == nView)
1156           iView = 0;
1157         if (iDet >= 0 && iDet < nDetWithZeros && iView >= 0 && iView < nView) {
1158           v[ix][iy] = ppcDetValue[iView][iDet].real();
1159           if (vImag)
1160             vImag[ix][iy] = ppcDetValue[iView][iDet].imag();
1161         } else {
1162           v[ix][iy] = 0;
1163           if (vImag)
1164             vImag[ix][iy] = 0;
1165         }
1166
1167       } else if (iInterpolationID == POLAR_INTERP_BILINEAR) {
1168         std::complex<double> vInterp = pBilinear->interpolate (ppdView[ix][iy], ppdDet[ix][iy]);
1169         v[ix][iy] = vInterp.real();
1170         if (vImag)
1171           vImag[ix][iy] = vInterp.imag();
1172       } else if (iInterpolationID == POLAR_INTERP_BICUBIC) {
1173         std::complex<double> vInterp = pBicubic->interpolate (ppdView[ix][iy], ppdDet[ix][iy]);
1174         v[ix][iy] = vInterp.real();
1175         if (vImag)
1176           vImag[ix][iy] = vInterp.imag();
1177       }
1178     }
1179   }
1180 }
1181
1182 bool
1183 Projections::initFromSomatomAR_STAR (int iNViews, int iNDets, unsigned char* pData, unsigned long lDataLength)
1184 {
1185   init (iNViews, iNDets);
1186   m_geometry = Scanner::GEOMETRY_EQUIANGULAR;
1187   m_dFocalLength = 510;
1188   m_dSourceDetectorLength = 890;
1189   m_detInc = convertDegreesToRadians (3.06976 / 60);
1190   m_dFanBeamAngle = iNDets * m_detInc;
1191   m_detStart = -(m_dFanBeamAngle / 2);
1192   m_rotInc = TWOPI / static_cast<double>(iNViews);
1193   m_rotStart = 0;
1194   m_dViewDiameter = sin (m_dFanBeamAngle / 2) * m_dFocalLength * 2;
1195
1196   if (! ((iNViews == 750 && lDataLength == 1560000L) || (iNViews == 950 && lDataLength == 1976000L) 
1197                 || (iNViews == 1500 && lDataLength == 3120000)))
1198     return false;
1199
1200   double dCenter = (iNDets - 1.) / 2.; // change from (Nm+1)/2 because of 0 vs. 1 indexing
1201   double* pdCosScale = new double [iNDets];
1202   for (int i = 0; i < iNDets; i++)
1203     pdCosScale[i] = 1. / cos ((i - dCenter) * m_detInc);
1204
1205   long lDataPos = 0;
1206   for (int iv = 0; iv < iNViews; iv++) {
1207     unsigned char* pArgBase = pData + lDataPos;
1208     unsigned char* p = pArgBase+0; SwapBytes4IfLittleEndian (p);
1209     // long lProjNumber = *reinterpret_cast<long*>(p);
1210
1211     p = pArgBase+20;  SwapBytes4IfLittleEndian (p);
1212     long lEscale = *reinterpret_cast<long*>(p);
1213
1214     p = pArgBase+28;  SwapBytes4IfLittleEndian (p);
1215     // long lTime = *reinterpret_cast<long*>(p);
1216
1217     p = pArgBase + 4; SwapBytes4IfLittleEndian (p);
1218     double dAlpha = *reinterpret_cast<float*>(p) + HALFPI;
1219
1220     p = pArgBase+12; SwapBytes4IfLittleEndian (p);
1221     // double dAlign = *reinterpret_cast<float*>(p);
1222
1223     p = pArgBase + 16; SwapBytes4IfLittleEndian (p);
1224     // double dMaxValue = *reinterpret_cast<float*>(p);
1225
1226     DetectorArray& detArray = getDetectorArray (iv);
1227     detArray.setViewAngle (dAlpha);
1228     DetectorValue* detval = detArray.detValues();
1229
1230     double dViewScale = 1. / (2294.4871 * ::pow (2.0, -lEscale));
1231     lDataPos += 32;
1232     for (int id = 0; id < iNDets; id++) {
1233       int iV = pData[lDataPos+1] + (pData[lDataPos] << 8);
1234       if (iV > 32767)   // two's complement signed conversion
1235         iV = iV - 65536;
1236       detval[id] = iV * dViewScale * pdCosScale[id];
1237       lDataPos += 2;
1238     }
1239 #if 1
1240     for (int k = iNDets - 2; k >= 0; k--)
1241       detval[k+1] = detval[k];
1242     detval[0] = 0;
1243 #endif
1244   }
1245
1246   delete pdCosScale;
1247   return true;
1248 }
1249
1250 Projections*
1251 Projections::interpolateToParallel () const
1252 {
1253   if (m_geometry == Scanner::GEOMETRY_PARALLEL)
1254     return const_cast<Projections*>(this);
1255
1256   int nDet = m_nDet;
1257   int nView = m_nView;
1258   Projections* pProjNew = new Projections (nView, nDet);
1259   pProjNew->m_geometry = Scanner::GEOMETRY_PARALLEL;
1260   pProjNew->m_dFocalLength = m_dFocalLength;
1261   pProjNew->m_dSourceDetectorLength = m_dSourceDetectorLength;
1262   pProjNew->m_dViewDiameter = m_dViewDiameter;
1263   pProjNew->m_dFanBeamAngle = m_dFanBeamAngle;
1264   pProjNew->m_calcTime  = 0;
1265   pProjNew->m_remark = m_remark;
1266   pProjNew->m_remark += "; Interpolate to Parallel";
1267   pProjNew->m_label.setLabelType (Array2dFileLabel::L_HISTORY);
1268   pProjNew->m_label.setLabelString (pProjNew->m_remark);
1269   pProjNew->m_label.setCalcTime (pProjNew->m_calcTime);
1270   pProjNew->m_label.setDateTime (pProjNew->m_year, pProjNew->m_month, pProjNew->m_day, pProjNew->m_hour, pProjNew->m_minute, pProjNew->m_second);
1271
1272   pProjNew->m_rotStart = 0;
1273 #ifdef CONVERT_PARALLEL_PI
1274   pProjNew->m_rotInc = PI / nView;;
1275 #else
1276   pProjNew->m_rotInc = TWOPI / nView;
1277 #endif
1278   pProjNew->m_detStart = -m_dViewDiameter / 2;
1279   pProjNew->m_detInc = m_dViewDiameter / nDet;
1280   if (isEven (nDet)) // even
1281     pProjNew->m_detInc = m_dViewDiameter / (nDet - 1);
1282
1283   ParallelRaysums parallel (this, ParallelRaysums::THETA_RANGE_NORMALIZE_TO_TWOPI);
1284
1285   double* pdThetaValuesForT = new double [pProjNew->nView()];
1286   double* pdRaysumsForT = new double [pProjNew->nView()];
1287
1288   // interpolate to evenly spaced theta (views)
1289   double dDetPos = pProjNew->m_detStart;
1290   for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc) {
1291       parallel.getThetaAndRaysumsForT (iD, pdThetaValuesForT, pdRaysumsForT);
1292
1293     double dViewAngle = m_rotStart;
1294     int iLastFloor = -1;
1295     LinearInterpolator<double> interp (pdThetaValuesForT, pdRaysumsForT, pProjNew->nView(), false);
1296     for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
1297       DetectorValue* detValues = pProjNew->getDetectorArray (iV).detValues();
1298       detValues[iD] = interp.interpolate (dViewAngle, &iLastFloor);
1299     }
1300   }
1301   delete pdThetaValuesForT;
1302   delete pdRaysumsForT;
1303
1304   // interpolate to evenly space t (detectors)
1305   double* pdOriginalDetPositions = new double [pProjNew->nDet()];
1306   parallel.getDetPositions (pdOriginalDetPositions);
1307
1308   double* pdDetValueCopy = new double [pProjNew->nDet()];
1309   double dViewAngle = m_rotStart;
1310   for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
1311     DetectorArray& detArray = pProjNew->getDetectorArray (iV);
1312     DetectorValue* detValues = detArray.detValues();
1313     detArray.setViewAngle (dViewAngle);
1314
1315     for (int i = 0; i < pProjNew->nDet(); i++)
1316       pdDetValueCopy[i] =   detValues[i];
1317
1318     double dDetPos = pProjNew->m_detStart;
1319     int iLastFloor = -1;
1320     LinearInterpolator<double> interp (pdOriginalDetPositions, pdDetValueCopy, pProjNew->nDet(), false);
1321     for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc)
1322       detValues[iD] = interp.interpolate (dDetPos, &iLastFloor);
1323   }
1324   delete pdDetValueCopy;
1325   delete pdOriginalDetPositions;
1326
1327   return pProjNew;
1328 }
1329
1330
1331 ///////////////////////////////////////////////////////////////////////////////
1332 //
1333 // Class ParallelRaysums
1334 //
1335 // Used for converting divergent beam raysums into Parallel raysums
1336 //
1337 ///////////////////////////////////////////////////////////////////////////////
1338
1339 ParallelRaysums::ParallelRaysums (const Projections* pProjections, int iThetaRange)
1340 : m_pCoordinates(NULL), m_iNumCoordinates(0), m_iNumView(pProjections->nView()), m_iNumDet(pProjections->nDet()),
1341   m_iThetaRange (iThetaRange)
1342 {
1343   int iGeometry = pProjections->geometry();
1344   double dDetInc = pProjections->detInc();
1345   double dDetStart = pProjections->detStart();
1346   double dFocalLength = pProjections->focalLength();
1347
1348   m_iNumCoordinates =  m_iNumView * m_iNumDet;
1349   m_pCoordinates = new ParallelRaysumCoordinate [m_iNumCoordinates];
1350   m_vecpCoordinates.reserve (m_iNumCoordinates);
1351   for (int i = 0; i < m_iNumCoordinates; i++)
1352     m_vecpCoordinates[i] = m_pCoordinates + i;
1353
1354   int iCoordinate = 0;
1355   for (int iV = 0; iV < m_iNumView; iV++) {
1356     double dViewAngle = pProjections->getDetectorArray(iV).viewAngle();
1357     const DetectorValue* detValues = pProjections->getDetectorArray(iV).detValues();
1358
1359     double dDetPos = dDetStart;
1360     for (int iD = 0; iD < m_iNumDet; iD++) {
1361       ParallelRaysumCoordinate* pC = m_vecpCoordinates[iCoordinate++];
1362
1363       if (iGeometry == Scanner::GEOMETRY_PARALLEL) {
1364         pC->m_dTheta = dViewAngle;
1365         pC->m_dT = dDetPos;
1366       } else if (iGeometry == Scanner::GEOMETRY_EQUILINEAR) {
1367         double dFanAngle = atan (dDetPos / pProjections->sourceDetectorLength());
1368         pC->m_dTheta = dViewAngle + dFanAngle;
1369         pC->m_dT = dFocalLength * sin(dFanAngle);        
1370
1371       } else if (iGeometry == Scanner::GEOMETRY_EQUIANGULAR) {
1372         // fan angle is same as dDetPos
1373         pC->m_dTheta = dViewAngle + dDetPos;
1374         pC->m_dT = dFocalLength * sin (dDetPos);        
1375       }
1376       if (m_iThetaRange != THETA_RANGE_UNCONSTRAINED) {
1377         pC->m_dTheta = normalizeAngle (pC->m_dTheta);
1378         if (m_iThetaRange == THETA_RANGE_FOLD_TO_PI && pC->m_dTheta >= PI) {
1379           pC->m_dTheta -= PI;
1380           pC->m_dT = -pC->m_dT;
1381         }
1382       }
1383       pC->m_dRaysum = detValues[iD];
1384       dDetPos += dDetInc;
1385     }
1386   }
1387 }
1388
1389 ParallelRaysums::~ParallelRaysums()
1390 {
1391   delete m_pCoordinates;
1392 }
1393
1394 ParallelRaysums::CoordinateContainer&
1395 ParallelRaysums::getSortedByTheta()
1396 {
1397   if (m_vecpSortedByTheta.size() == 0) {
1398     m_vecpSortedByTheta.resize (m_iNumCoordinates);
1399     for (int i = 0; i < m_iNumCoordinates; i++)
1400       m_vecpSortedByTheta[i] = m_vecpCoordinates[i];
1401     std::sort (m_vecpSortedByTheta.begin(), m_vecpSortedByTheta.end(), ParallelRaysumCoordinate::compareByTheta);
1402   }
1403
1404   return m_vecpSortedByTheta;
1405 }
1406
1407 ParallelRaysums::CoordinateContainer&
1408 ParallelRaysums::getSortedByT()
1409 {
1410   if (m_vecpSortedByT.size() == 0) {
1411     m_vecpSortedByT.resize (m_iNumCoordinates);
1412     for (int i = 0; i < m_iNumCoordinates; i++)
1413       m_vecpSortedByT[i] = m_vecpCoordinates[i];
1414     std::sort (m_vecpSortedByT.begin(), m_vecpSortedByT.end(), ParallelRaysumCoordinate::compareByT);
1415   }
1416
1417   return m_vecpSortedByT;
1418 }
1419
1420
1421 void
1422 ParallelRaysums::getLimits (double* dMinT, double* dMaxT, double* dMinTheta, double* dMaxTheta) const
1423 {
1424   if (m_iNumCoordinates <= 0)
1425     return;
1426
1427   *dMinT = *dMaxT = m_vecpCoordinates[0]->m_dT;
1428   *dMinTheta = *dMaxTheta = m_vecpCoordinates[0]->m_dTheta;
1429
1430   for (int i = 0; i < m_iNumCoordinates; i++) {
1431     double dT = m_vecpCoordinates[i]->m_dT;
1432     double dTheta = m_vecpCoordinates[i]->m_dTheta;
1433
1434     if (dT < *dMinT)
1435       *dMinT = dT;
1436     else if (dT > *dMaxT)
1437       *dMaxT = dT;
1438
1439     if (dTheta < *dMinTheta)
1440       *dMinTheta = dTheta;
1441     else if (dTheta > *dMaxTheta)
1442       *dMaxTheta = dTheta;
1443   }
1444 }
1445
1446 void
1447 ParallelRaysums::getThetaAndRaysumsForT (int iTheta, double* pTheta, double* pRaysum)
1448 {
1449   const CoordinateContainer& coordsT = getSortedByT();
1450
1451   int iBase = iTheta * m_iNumView;
1452   for (int i = 0; i < m_iNumView; i++) {
1453     int iPos = iBase + i;
1454     pTheta[i] = coordsT[iPos]->m_dTheta;
1455     pRaysum[i] = coordsT[iPos]->m_dRaysum;
1456   }
1457 }
1458
1459 void
1460 ParallelRaysums::getDetPositions (double* pdDetPos)
1461 {
1462   const CoordinateContainer& coordsT = getSortedByT();
1463
1464   int iPos = 0;
1465   for (int i = 0; i < m_iNumDet; i++) {
1466     pdDetPos[i] = coordsT[iPos]->m_dT;
1467     iPos += m_iNumView;
1468   }
1469 }