1 /*****************************************************************************
4 ** Name: projections.cpp Projection data classes
5 ** Programmer: Kevin Rosenberg
6 ** Date Started: Aug 84
8 ** This is part of the CTSim program
9 ** Copyright (c) 1983-2001 Kevin Rosenberg
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.
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.
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 ******************************************************************************/
29 #include "interpolator.h"
31 const kuint16 Projections::m_signature = ('P'*256 + 'J');
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;
38 const char* const Projections::s_aszInterpName[] =
45 const char* const Projections::s_aszInterpTitle[] =
52 const int Projections::s_iInterpCount = sizeof(s_aszInterpName) / sizeof(char*);
57 * Projections Constructor for projections matrix storage
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
67 Projections::Projections (const Scanner& scanner)
70 initFromScanner (scanner);
74 Projections::Projections (const int nView, const int nDet)
80 Projections::Projections (void)
86 Projections::~Projections (void)
92 Projections::convertInterpNameToID (const char* const interpName)
94 int interpID = POLAR_INTERP_INVALID;
96 for (int i = 0; i < s_iInterpCount; i++)
97 if (strcasecmp (interpName, s_aszInterpName[i]) == 0) {
106 Projections::convertInterpIDToName (const int interpID)
108 static const char *interpName = "";
110 if (interpID >= 0 && interpID < s_iInterpCount)
111 return (s_aszInterpName[interpID]);
117 Projections::convertInterpIDToTitle (const int interpID)
119 static const char *interpTitle = "";
121 if (interpID >= 0 && interpID < s_iInterpCount)
122 return (s_aszInterpTitle[interpID]);
124 return (interpTitle);
130 Projections::init (const int nView, const int nDet)
132 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
137 time_t t = time (NULL);
138 tm* lt = localtime (&t);
139 m_year = lt->tm_year;
140 m_month = lt->tm_mon;
142 m_hour = lt->tm_hour;
143 m_minute = lt->tm_min;
144 m_second = lt->tm_sec;
148 Projections::initFromScanner (const Scanner& scanner)
150 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
152 init (scanner.nView(), scanner.nDet());
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();
166 Projections::setNView (int nView) // used by MPI to reduce # of views
169 init (nView, m_nDet);
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.
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
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.
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,
203 Projections::Helical180LI(int interpolation_view)
205 if (m_geometry == Scanner::GEOMETRY_INVALID)
207 std::cerr << "Invalid geometry " << m_geometry << std::endl;
210 else if (m_geometry == Scanner::GEOMETRY_PARALLEL)
212 std::cerr << "Helical 180LI not yet implemented for PARALLEL geometry"
216 else if (m_geometry == Scanner::GEOMETRY_EQUILINEAR)
218 std::cerr << "Helical 180LI not yet implemented for EQUILINEAR geometry"
222 else if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR)
224 return Helical180LI_Equiangular(interpolation_view);
228 std::cerr << "Invalid geometry in projection data file" << m_geometry
234 Projections::Helical180LI_Equiangular(int interpView)
236 double dbeta = m_rotInc;
237 double dgamma = m_detInc;
238 double fanAngle = m_dFanBeamAngle;
241 // is there enough data in the data set? Should have 2(Pi+fanAngle)
243 if ( m_nView < static_cast<int>((2*( PI + fanAngle ) ) / dbeta) -1 ){
244 std::cerr << "Data set does not include 360 +2*FanBeamAngle views"
249 if (interpView < 0) // use default position at PI+fanAngle
251 interpView = static_cast<int> ((PI+fanAngle)/dbeta);
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)
260 std::cerr << "There isn't PI+fanAngle of data on either side of the requested interpolation view" << std::endl;
263 offsetView = interpView - static_cast<int>((PI+fanAngle)/dbeta);
266 int last_interp_view = static_cast<int> ((PI+fanAngle)/dbeta);
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++)
280 int last_acq_view = 2*last_interp_view;
281 for ( int iView = 0 ; iView <= last_acq_view; iView++) {
282 double beta = iView * dbeta;
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 )
293 else // (beta > PI+fanAngle)
295 //newbeta = beta +2*gamma - 180;
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);
304 //std::cout << beta << " "<< gamma << " " << newbeta << " " << newgamma <<" " << iView-offsetView << " " << iDet << " " << newiView << " " << newiDet << std::endl;
305 //std::cout << iView-offsetView << " " << iDet << " " << newiView << " " << newiDet << std::endl;
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)
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)
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)
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)
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)
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)
350 ; // outside region of interest
356 m_projData = newdetarray;
357 m_nView = last_interp_view+1;
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.
372 Projections::HalfScanFeather(void)
374 double dbeta = m_rotInc;
375 double dgamma = m_detInc;
376 double fanAngle = m_dFanBeamAngle;
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;
383 if (m_geometry == Scanner::GEOMETRY_INVALID) {
384 std::cerr << "Invalid geometry " << m_geometry << std::endl;
388 if (m_geometry == Scanner::GEOMETRY_PARALLEL) {
389 std::cerr << "HalfScanFeather not yet implemented for PARALLEL geometry"<< std::endl;
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
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);
405 DetectorValue* detval2 = (m_projData[iView2])->detValues();
406 DetectorValue* detval1 = (m_projData[iView1])->detValues();
408 detval1[iDet1] = detval2[iDet2] ;
410 double x, w1,w2,beta1, gamma1;
413 if ( beta1 <= (fanAngle - 2*gamma1) )
414 x = beta1 / ( fanAngle - 2*gamma1);
415 else if ( (fanAngle - 2*gamma1 <= beta1 ) && beta1 <= PI - 2*gamma1)
417 else if ( (PI - 2*gamma1 <= beta1 ) && ( beta1 <=PI + fanAngle) )
418 x = (PI +fanAngle - beta1)/(fanAngle + 2*gamma1);
420 std::cerr << "Shouldn't be here!"<< std::endl;
423 w1 = (3*x - 2*x*x)*x;
425 detval1[iDet1] *= w1;
426 detval2[iDet2] *= w2;
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;
447 Projections::newProjData (void)
450 sys_error(ERR_WARNING, "m_projData != NULL [newProjData]");
452 if (m_nView > 0 && m_nDet) {
453 m_projData = new DetectorArray* [m_nView];
455 for (int i = 0; i < m_nView; i++)
456 m_projData[i] = new DetectorArray (m_nDet);
462 * projections_free Free memory allocated to projections
465 * projections_free(proj)
466 * Projections& proj Projectionss to be deallocated
470 Projections::deleteProjData (void)
472 if (m_projData != NULL) {
473 for (int i = 0; i < m_nView; i++)
474 delete m_projData[i];
483 * Projections::headerWwrite Write data header for projections file
488 Projections::headerWrite (fnetorderstream& fs)
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;
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;
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);
540 m_headerSize = fs.tellp();
541 _hsize = m_headerSize;
543 fs.writeInt16 (_hsize);
551 * projections_read_header Read data header for projections file
555 Projections::headerRead (fnetorderstream& fs)
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;
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);
582 fs.readInt16 (_hour);
583 fs.readInt16 (_minute);
584 fs.readInt16 (_second);
585 fs.readInt16 (_remarksize);
588 sys_error (ERR_SEVERE, "Error reading header information , _remarksize=%d [projections_read_header]", _remarksize);
592 if (_signature != m_signature) {
593 sys_error (ERR_SEVERE, "File %s does not have a valid projection file signature", m_filename.c_str());
597 char* pszRemarkStorage = new char [_remarksize+1];
598 fs.read (pszRemarkStorage, _remarksize);
600 sys_error (ERR_SEVERE, "Error reading remark, _remarksize = %d", _remarksize);
603 pszRemarkStorage[_remarksize] = 0;
604 m_remark = pszRemarkStorage;
605 delete pszRemarkStorage;
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);
613 m_headerSize = _hsize;
617 m_calcTime = _calcTime;
618 m_rotStart = _rotStart;
620 m_detStart = _detStart;
622 m_dFocalLength = _focalLength;
623 m_dSourceDetectorLength = _sourceDetectorLength;
624 m_dViewDiameter = _viewDiameter;
625 m_dFanBeamAngle = _fanBeamAngle;
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);
642 Projections::read (const std::string& filename)
644 return read (filename.c_str());
649 Projections::read (const char* filename)
651 m_filename = filename;
653 frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary);
655 frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary); // | std::ios::nocreate);
661 if (! headerRead (fileRead))
667 for (int i = 0; i < m_nView; i++) {
668 if (! detarrayRead (fileRead, *m_projData[i], i))
678 Projections::copyViewData (const std::string& filename, std::ostream& os, int startView, int endView)
680 return copyViewData (filename.c_str(), os, startView, endView);
684 Projections::copyViewData (const char* const filename, std::ostream& os, int startView, int endView)
686 frnetorderstream is (filename, std::ios::in | std::ios::binary);
687 kuint16 sizeHeader, signature;
688 kuint32 _nView, _nDet;
692 sys_error (ERR_SEVERE, "Unable to read projection file %s", filename);
696 is.readInt16 (sizeHeader);
697 is.readInt16 (signature);
698 is.readInt32 (_nView);
699 is.readInt32 (_nDet);
703 if (signature != m_signature) {
704 sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
710 if (startView > nView - 1)
712 if (endView < 0 || endView > nView - 1)
715 if (startView > endView) { // swap if start > end
716 int tempView = endView;
718 startView = tempView;
721 int sizeView = 8 /* view_angle */ + 4 /* nDet */ + (4 * nDet);
722 unsigned char* pViewData = new unsigned char [sizeView];
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())
734 sys_error (ERR_SEVERE, "Error reading projection file");
736 sys_error (ERR_SEVERE, "Error writing projection file");
738 return (! (is.fail() | os.fail()));
742 Projections::copyHeader (const std::string& filename, std::ostream& os)
744 return copyHeader (filename.c_str(), os);
748 Projections::copyHeader (const char* const filename, std::ostream& os)
750 frnetorderstream is (filename, std::ios::in | std::ios::binary);
751 kuint16 sizeHeader, signature;
752 is.readInt16 (sizeHeader);
753 is.readInt16 (signature);
755 if (signature != m_signature) {
756 sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
760 unsigned char* pHdrData = new unsigned char [sizeHeader];
761 is.read (reinterpret_cast<char*>(pHdrData), sizeHeader);
763 sys_error (ERR_SEVERE, "Error reading header");
767 os.write (reinterpret_cast<char*>(pHdrData), sizeHeader);
769 sys_error (ERR_SEVERE, "Error writing header");
777 Projections::write (const std::string& filename)
779 return write (filename.c_str());
783 Projections::write (const char* filename)
785 frnetorderstream fs (filename, std::ios::out | std::ios::binary | std::ios::trunc | std::ios::ate);
786 m_filename = filename;
788 sys_error (ERR_SEVERE, "Error opening file %s for output [projections_create]", filename);
792 if (! headerWrite (fs))
795 if (m_projData != NULL) {
796 for (int i = 0; i < m_nView; i++) {
797 if (! detarrayWrite (fs, *m_projData[i], i))
810 * detarrayRead Read a Detector Array structure from the disk
813 * detarrayRead (proj, darray, view_num)
814 * DETARRAY *darray Detector array storage location to be filled
815 * int view_num View number to read
819 Projections::detarrayRead (fnetorderstream& fs, DetectorArray& darray, const int iview)
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();
829 fs.seekg (start_data);
831 fs.readFloat64 (view_angle);
833 darray.setViewAngle (view_angle);
834 // darray.setNDet ( nDet);
836 for (unsigned int i = 0; i < nDet; i++) {
838 fs.readFloat32 (detval);
839 detval_ptr[i] = detval;
849 * detarrayWrite Write detector array data to the disk
852 * detarrayWrite (darray, view_num)
853 * DETARRAY *darray Detector array data to be written
854 * int view_num View number to write
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
863 Projections::detarrayWrite (fnetorderstream& fs, const DetectorArray& darray, const int iview)
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();
873 fs.seekp (start_data);
875 sys_error (ERR_SEVERE, "Error seeking detectory array [detarrayWrite]");
879 fs.writeFloat64 (view_angle);
880 fs.writeInt32 (nDet);
882 for (unsigned int i = 0; i < nDet; i++) {
883 kfloat32 detval = detval_ptr[i];
884 fs.writeFloat32 (detval);
894 * printProjectionData Print projections data
897 * printProjectionData ()
901 Projections::printProjectionData ()
903 printProjectionData (0, nView() - 1);
907 Projections::printProjectionData (int startView, int endView)
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) {
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]);
937 Projections::printScanInfo (std::ostringstream& os) const
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";
955 Projections::convertPolar (ImageFile& rIF, int iInterpolationID)
957 unsigned int nx = rIF.nx();
958 unsigned int ny = rIF.ny();
959 ImageFileArray v = rIF.getArray();
960 ImageFileArray vImag = rIF.getImaginaryArray();
962 if (! v || nx == 0 || ny == 0)
965 Projections* pProj = this;
966 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
967 pProj = interpolateToParallel();
969 Array2d<double> adView (nx, ny);
970 Array2d<double> adDet (nx, ny);
971 double** ppdView = adView.getArray();
972 double** ppdDet = adDet.getArray();
974 std::complex<double>** ppcDetValue = new std::complex<double>* [pProj->m_nView];
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);
983 pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet, pProj->m_nDet, 1., pProj->m_detInc);
985 pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet,
986 pProj->m_nDet, iInterpolationID);
988 for (iView = 0; iView < pProj->m_nView; iView++)
989 delete [] ppcDetValue[iView];
990 delete [] ppcDetValue;
992 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
1000 Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad)
1003 rIF.arrayDataClear();
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();
1013 if (! v || nx == 0 || ny == 0)
1016 Projections* pProj = this;
1017 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
1018 pProj = interpolateToParallel();
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.05);
1023 double dZeropadRatio = static_cast<double>(iNumInterpDetWithZeros) / static_cast<double>(iInterpDet);
1025 fftw_complex* pcIn = static_cast<fftw_complex*> (fftw_malloc (sizeof(fftw_complex) * iNumInterpDetWithZeros));
1026 fftw_plan plan = fftw_plan_dft_1d (iNumInterpDetWithZeros, pcIn, pcIn, FFTW_FORWARD, FFTW_ESTIMATE);
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);
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;
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][0] = projInterp.interpolate (dInterpPos) * dProjScale;
1047 Fourier::shuffleFourierToNaturalOrder (pcIn, iInterpDet);
1048 if (iZerosAdded > 0) {
1049 for (int iDet1 = iInterpDet -1; iDet1 >= iMidPoint; iDet1--) {
1050 pcIn[iDet1+iZerosAdded][0] = pcIn[iDet1][0];
1051 pcIn[iDet1+iZerosAdded][1] = pcIn[iDet1][1];
1053 for (int iDet2 = iMidPoint; iDet2 < iMidPoint + iZerosAdded; iDet2++)
1054 pcIn[iDet2][0] = pcIn[iDet2][1] = 0;
1057 fftw_execute (plan);
1059 ppcDetValue[iView] = new std::complex<double> [iNumInterpDetWithZeros];
1060 for (int iD = 0; iD < iNumInterpDetWithZeros; iD++) {
1061 ppcDetValue[iView][iD] = std::complex<double> (pcIn[iD][0] * dFFTScale, pcIn[iD][1] * dFFTScale);
1064 Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], iNumInterpDetWithZeros);
1068 fftw_destroy_plan (plan);
1070 Array2d<double> adView (nx, ny);
1071 Array2d<double> adDet (nx, ny);
1072 double** ppdView = adView.getArray();
1073 double** ppdDet = adDet.getArray();
1074 pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet, iNumInterpDetWithZeros, dZeropadRatio,
1075 pProj->m_detInc * dInterpScale);
1077 pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet,
1078 iNumInterpDetWithZeros, iInterpolationID);
1080 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
1083 for (int i = 0; i < m_nView; i++)
1084 delete [] ppcDetValue[i];
1085 delete [] ppcDetValue;
1093 Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double** ppdView, double** ppdDet,
1094 int iNumDetWithZeros, double dZeropadRatio, double dDetInc)
1096 double dLength = viewDiameter();
1097 double xMin = -dLength / 2;
1098 double xMax = xMin + dLength;
1099 double yMin = -dLength / 2;
1100 double yMax = yMin + dLength;
1101 double xCent = (xMin + xMax) / 2;
1102 double yCent = (yMin + yMax) / 2;
1104 xMin = (xMin - xCent) * dZeropadRatio + xCent;
1105 xMax = (xMax - xCent) * dZeropadRatio + xCent;
1106 yMin = (yMin - yCent) * dZeropadRatio + yCent;
1107 yMax = (yMax - yCent) * dZeropadRatio + yCent;
1109 double xInc = (xMax - xMin) / nx; // size of cells
1110 double yInc = (yMax - yMin) / ny;
1112 double dDetCenter = (iNumDetWithZeros - 1) / 2.; // index refering to L=0 projection
1113 // +1 is correct for frequency data, ndet-1 is correct for projections
1114 // if (isEven (iNumDetWithZeros))
1115 // dDetCenter = (iNumDetWithZeros + 0) / 2;
1117 // Calculates polar coordinates (view#, det#) for each point on phantom grid
1118 double x = xMin + xInc / 2; // Rectang coords of center of pixel
1119 for (unsigned int ix = 0; ix < nx; x += xInc, ix++) {
1120 double y = yMin + yInc / 2;
1121 for (unsigned int iy = 0; iy < ny; y += yInc, iy++) {
1122 double r = ::sqrt (x * x + y * y);
1123 double phi = atan2 (y, x);
1125 if (phi <= -m_rotInc / 2)
1127 if (phi >= PI - (m_rotInc / 2)) {
1132 ppdView[ix][iy] = (phi - m_rotStart) / m_rotInc;
1133 ppdDet[ix][iy] = (r / dDetInc) + dDetCenter;
1139 Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag,
1140 unsigned int nx, unsigned int ny, std::complex<double>** ppcDetValue, double** ppdView,
1141 double** ppdDet, unsigned int nView, unsigned int nDet, unsigned int nDetWithZeros, int iInterpolationID)
1143 typedef std::complex<double> complexValue;
1145 BilinearPolarInterpolator<complexValue>* pBilinear = NULL;
1146 BicubicPolyInterpolator<complexValue>* pBicubic = NULL;
1147 if (iInterpolationID == POLAR_INTERP_BILINEAR)
1148 pBilinear = new BilinearPolarInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
1149 else if (iInterpolationID == POLAR_INTERP_BICUBIC)
1150 pBicubic = new BicubicPolyInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
1152 for (unsigned int ix = 0; ix < ny; ix++) {
1153 for (unsigned int iy = 0; iy < ny; iy++) {
1154 if (iInterpolationID == POLAR_INTERP_NEAREST) {
1155 unsigned int iView = nearest<int> (ppdView[ix][iy]);
1156 unsigned int iDet = nearest<int> (ppdDet[ix][iy]);
1159 if (iDet >= 0 && iDet < nDetWithZeros && iView >= 0 && iView < nView) {
1160 v[ix][iy] = ppcDetValue[iView][iDet].real();
1162 vImag[ix][iy] = ppcDetValue[iView][iDet].imag();
1169 } else if (iInterpolationID == POLAR_INTERP_BILINEAR) {
1170 std::complex<double> vInterp = pBilinear->interpolate (ppdView[ix][iy], ppdDet[ix][iy]);
1171 v[ix][iy] = vInterp.real();
1173 vImag[ix][iy] = vInterp.imag();
1174 } else if (iInterpolationID == POLAR_INTERP_BICUBIC) {
1175 std::complex<double> vInterp = pBicubic->interpolate (ppdView[ix][iy], ppdDet[ix][iy]);
1176 v[ix][iy] = vInterp.real();
1178 vImag[ix][iy] = vInterp.imag();
1185 Projections::initFromSomatomAR_STAR (int iNViews, int iNDets, unsigned char* pData, unsigned long lDataLength)
1187 init (iNViews, iNDets);
1188 m_geometry = Scanner::GEOMETRY_EQUIANGULAR;
1189 m_dFocalLength = 510;
1190 m_dSourceDetectorLength = 890;
1191 m_detInc = convertDegreesToRadians (3.06976 / 60);
1192 m_dFanBeamAngle = iNDets * m_detInc;
1193 m_detStart = -(m_dFanBeamAngle / 2);
1194 m_rotInc = TWOPI / static_cast<double>(iNViews);
1196 m_dViewDiameter = sin (m_dFanBeamAngle / 2) * m_dFocalLength * 2;
1198 if (! ((iNViews == 750 && lDataLength == 1560000L) || (iNViews == 950 && lDataLength == 1976000L)
1199 || (iNViews == 1500 && lDataLength == 3120000)))
1202 double dCenter = (iNDets - 1.) / 2.; // change from (Nm+1)/2 because of 0 vs. 1 indexing
1203 double* pdCosScale = new double [iNDets];
1204 for (int i = 0; i < iNDets; i++)
1205 pdCosScale[i] = 1. / cos ((i - dCenter) * m_detInc);
1208 for (int iv = 0; iv < iNViews; iv++) {
1209 unsigned char* pArgBase = pData + lDataPos;
1210 unsigned char* p = pArgBase+0; SwapBytes4IfLittleEndian (p);
1211 // long lProjNumber = *reinterpret_cast<long*>(p);
1213 p = pArgBase+20; SwapBytes4IfLittleEndian (p);
1214 long lEscale = *reinterpret_cast<long*>(p);
1216 p = pArgBase+28; SwapBytes4IfLittleEndian (p);
1217 // long lTime = *reinterpret_cast<long*>(p);
1219 p = pArgBase + 4; SwapBytes4IfLittleEndian (p);
1220 double dAlpha = *reinterpret_cast<float*>(p) + HALFPI;
1222 p = pArgBase+12; SwapBytes4IfLittleEndian (p);
1223 // double dAlign = *reinterpret_cast<float*>(p);
1225 p = pArgBase + 16; SwapBytes4IfLittleEndian (p);
1226 // double dMaxValue = *reinterpret_cast<float*>(p);
1228 DetectorArray& detArray = getDetectorArray (iv);
1229 detArray.setViewAngle (dAlpha);
1230 DetectorValue* detval = detArray.detValues();
1232 double dViewScale = 1. / (2294.4871 * ::pow (2.0, -lEscale));
1234 for (int id = 0; id < iNDets; id++) {
1235 int iV = pData[lDataPos+1] + (pData[lDataPos] << 8);
1236 if (iV > 32767) // two's complement signed conversion
1238 detval[id] = iV * dViewScale * pdCosScale[id];
1242 for (int k = iNDets - 2; k >= 0; k--)
1243 detval[k+1] = detval[k];
1253 Projections::interpolateToParallel () const
1255 if (m_geometry == Scanner::GEOMETRY_PARALLEL)
1256 return const_cast<Projections*>(this);
1259 int nView = m_nView;
1260 Projections* pProjNew = new Projections (nView, nDet);
1261 pProjNew->m_geometry = Scanner::GEOMETRY_PARALLEL;
1262 pProjNew->m_dFocalLength = m_dFocalLength;
1263 pProjNew->m_dSourceDetectorLength = m_dSourceDetectorLength;
1264 pProjNew->m_dViewDiameter = m_dViewDiameter;
1265 pProjNew->m_dFanBeamAngle = m_dFanBeamAngle;
1266 pProjNew->m_calcTime = 0;
1267 pProjNew->m_remark = m_remark;
1268 pProjNew->m_remark += "; Interpolate to Parallel";
1269 pProjNew->m_label.setLabelType (Array2dFileLabel::L_HISTORY);
1270 pProjNew->m_label.setLabelString (pProjNew->m_remark);
1271 pProjNew->m_label.setCalcTime (pProjNew->m_calcTime);
1272 pProjNew->m_label.setDateTime (pProjNew->m_year, pProjNew->m_month, pProjNew->m_day, pProjNew->m_hour, pProjNew->m_minute, pProjNew->m_second);
1274 pProjNew->m_rotStart = 0;
1275 #ifdef CONVERT_PARALLEL_PI
1276 pProjNew->m_rotInc = PI / nView;;
1278 pProjNew->m_rotInc = TWOPI / nView;
1280 pProjNew->m_detStart = -m_dViewDiameter / 2;
1281 pProjNew->m_detInc = m_dViewDiameter / nDet;
1282 if (isEven (nDet)) // even
1283 pProjNew->m_detInc = m_dViewDiameter / (nDet - 1);
1285 ParallelRaysums parallel (this, ParallelRaysums::THETA_RANGE_NORMALIZE_TO_TWOPI);
1287 double* pdThetaValuesForT = new double [pProjNew->nView()];
1288 double* pdRaysumsForT = new double [pProjNew->nView()];
1290 // interpolate to evenly spaced theta (views)
1291 double dDetPos = pProjNew->m_detStart;
1292 for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc) {
1293 parallel.getThetaAndRaysumsForT (iD, pdThetaValuesForT, pdRaysumsForT);
1295 double dViewAngle = m_rotStart;
1296 int iLastFloor = -1;
1297 LinearInterpolator<double> interp (pdThetaValuesForT, pdRaysumsForT, pProjNew->nView(), false);
1298 for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
1299 DetectorValue* detValues = pProjNew->getDetectorArray (iV).detValues();
1300 detValues[iD] = interp.interpolate (dViewAngle, &iLastFloor);
1303 delete pdThetaValuesForT;
1304 delete pdRaysumsForT;
1306 // interpolate to evenly space t (detectors)
1307 double* pdOriginalDetPositions = new double [pProjNew->nDet()];
1308 parallel.getDetPositions (pdOriginalDetPositions);
1310 double* pdDetValueCopy = new double [pProjNew->nDet()];
1311 double dViewAngle = m_rotStart;
1312 for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
1313 DetectorArray& detArray = pProjNew->getDetectorArray (iV);
1314 DetectorValue* detValues = detArray.detValues();
1315 detArray.setViewAngle (dViewAngle);
1317 for (int i = 0; i < pProjNew->nDet(); i++)
1318 pdDetValueCopy[i] = detValues[i];
1320 double dDetPos = pProjNew->m_detStart;
1321 int iLastFloor = -1;
1322 LinearInterpolator<double> interp (pdOriginalDetPositions, pdDetValueCopy, pProjNew->nDet(), false);
1323 for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc)
1324 detValues[iD] = interp.interpolate (dDetPos, &iLastFloor);
1326 delete pdDetValueCopy;
1327 delete pdOriginalDetPositions;
1333 ///////////////////////////////////////////////////////////////////////////////
1335 // Class ParallelRaysums
1337 // Used for converting divergent beam raysums into Parallel raysums
1339 ///////////////////////////////////////////////////////////////////////////////
1341 ParallelRaysums::ParallelRaysums (const Projections* pProjections, int iThetaRange)
1342 : m_pCoordinates(NULL), m_iNumCoordinates(0), m_iNumView(pProjections->nView()), m_iNumDet(pProjections->nDet()),
1343 m_iThetaRange (iThetaRange)
1345 int iGeometry = pProjections->geometry();
1346 double dDetInc = pProjections->detInc();
1347 double dDetStart = pProjections->detStart();
1348 double dFocalLength = pProjections->focalLength();
1350 m_iNumCoordinates = m_iNumView * m_iNumDet;
1351 m_pCoordinates = new ParallelRaysumCoordinate [m_iNumCoordinates];
1352 m_vecpCoordinates.reserve (m_iNumCoordinates);
1353 for (int i = 0; i < m_iNumCoordinates; i++)
1354 m_vecpCoordinates[i] = m_pCoordinates + i;
1356 int iCoordinate = 0;
1357 for (int iV = 0; iV < m_iNumView; iV++) {
1358 double dViewAngle = pProjections->getDetectorArray(iV).viewAngle();
1359 const DetectorValue* detValues = pProjections->getDetectorArray(iV).detValues();
1361 double dDetPos = dDetStart;
1362 for (int iD = 0; iD < m_iNumDet; iD++) {
1363 ParallelRaysumCoordinate* pC = m_vecpCoordinates[iCoordinate++];
1365 if (iGeometry == Scanner::GEOMETRY_PARALLEL) {
1366 pC->m_dTheta = dViewAngle;
1368 } else if (iGeometry == Scanner::GEOMETRY_EQUILINEAR) {
1369 double dFanAngle = atan (dDetPos / pProjections->sourceDetectorLength());
1370 pC->m_dTheta = dViewAngle + dFanAngle;
1371 pC->m_dT = dFocalLength * sin(dFanAngle);
1373 } else if (iGeometry == Scanner::GEOMETRY_EQUIANGULAR) {
1374 // fan angle is same as dDetPos
1375 pC->m_dTheta = dViewAngle + dDetPos;
1376 pC->m_dT = dFocalLength * sin (dDetPos);
1378 if (m_iThetaRange != THETA_RANGE_UNCONSTRAINED) {
1379 pC->m_dTheta = normalizeAngle (pC->m_dTheta);
1380 if (m_iThetaRange == THETA_RANGE_FOLD_TO_PI && pC->m_dTheta >= PI) {
1382 pC->m_dT = -pC->m_dT;
1385 pC->m_dRaysum = detValues[iD];
1391 ParallelRaysums::~ParallelRaysums()
1393 delete m_pCoordinates;
1396 ParallelRaysums::CoordinateContainer&
1397 ParallelRaysums::getSortedByTheta()
1399 if (m_vecpSortedByTheta.size() == 0) {
1400 m_vecpSortedByTheta.resize (m_iNumCoordinates);
1401 for (int i = 0; i < m_iNumCoordinates; i++)
1402 m_vecpSortedByTheta[i] = m_vecpCoordinates[i];
1403 std::sort (m_vecpSortedByTheta.begin(), m_vecpSortedByTheta.end(), ParallelRaysumCoordinate::compareByTheta);
1406 return m_vecpSortedByTheta;
1409 ParallelRaysums::CoordinateContainer&
1410 ParallelRaysums::getSortedByT()
1412 if (m_vecpSortedByT.size() == 0) {
1413 m_vecpSortedByT.resize (m_iNumCoordinates);
1414 for (int i = 0; i < m_iNumCoordinates; i++)
1415 m_vecpSortedByT[i] = m_vecpCoordinates[i];
1416 std::sort (m_vecpSortedByT.begin(), m_vecpSortedByT.end(), ParallelRaysumCoordinate::compareByT);
1419 return m_vecpSortedByT;
1424 ParallelRaysums::getLimits (double* dMinT, double* dMaxT, double* dMinTheta, double* dMaxTheta) const
1426 if (m_iNumCoordinates <= 0)
1429 *dMinT = *dMaxT = m_vecpCoordinates[0]->m_dT;
1430 *dMinTheta = *dMaxTheta = m_vecpCoordinates[0]->m_dTheta;
1432 for (int i = 0; i < m_iNumCoordinates; i++) {
1433 double dT = m_vecpCoordinates[i]->m_dT;
1434 double dTheta = m_vecpCoordinates[i]->m_dTheta;
1438 else if (dT > *dMaxT)
1441 if (dTheta < *dMinTheta)
1442 *dMinTheta = dTheta;
1443 else if (dTheta > *dMaxTheta)
1444 *dMaxTheta = dTheta;
1449 ParallelRaysums::getThetaAndRaysumsForT (int iTheta, double* pTheta, double* pRaysum)
1451 const CoordinateContainer& coordsT = getSortedByT();
1453 int iBase = iTheta * m_iNumView;
1454 for (int i = 0; i < m_iNumView; i++) {
1455 int iPos = iBase + i;
1456 pTheta[i] = coordsT[iPos]->m_dTheta;
1457 pRaysum[i] = coordsT[iPos]->m_dRaysum;
1462 ParallelRaysums::getDetPositions (double* pdDetPos)
1464 const CoordinateContainer& coordsT = getSortedByT();
1467 for (int i = 0; i < m_iNumDet; i++) {
1468 pdDetPos[i] = coordsT[iPos]->m_dT;