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