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
11 ** $Id: projections.cpp,v 1.80 2002/06/27 03:19:23 kevin Exp $
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 ******************************************************************************/
30 const kuint16 Projections::m_signature = ('P'*256 + 'J');
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;
37 const char* const Projections::s_aszInterpName[] =
44 const char* const Projections::s_aszInterpTitle[] =
51 const int Projections::s_iInterpCount = sizeof(s_aszInterpName) / sizeof(char*);
56 * Projections Constructor for projections matrix storage
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
66 Projections::Projections (const Scanner& scanner)
69 initFromScanner (scanner);
73 Projections::Projections (const int nView, const int nDet)
79 Projections::Projections (void)
85 Projections::~Projections (void)
91 Projections::convertInterpNameToID (const char* const interpName)
93 int interpID = POLAR_INTERP_INVALID;
95 for (int i = 0; i < s_iInterpCount; i++)
96 if (strcasecmp (interpName, s_aszInterpName[i]) == 0) {
105 Projections::convertInterpIDToName (const int interpID)
107 static const char *interpName = "";
109 if (interpID >= 0 && interpID < s_iInterpCount)
110 return (s_aszInterpName[interpID]);
116 Projections::convertInterpIDToTitle (const int interpID)
118 static const char *interpTitle = "";
120 if (interpID >= 0 && interpID < s_iInterpCount)
121 return (s_aszInterpTitle[interpID]);
123 return (interpTitle);
129 Projections::init (const int nView, const int nDet)
131 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
136 time_t t = time (NULL);
137 tm* lt = localtime (&t);
138 m_year = lt->tm_year;
139 m_month = lt->tm_mon;
141 m_hour = lt->tm_hour;
142 m_minute = lt->tm_min;
143 m_second = lt->tm_sec;
147 Projections::initFromScanner (const Scanner& scanner)
149 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
151 init (scanner.nView(), scanner.nDet());
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();
165 Projections::setNView (int nView) // used by MPI to reduce # of views
168 init (nView, m_nDet);
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.
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
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.
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,
202 Projections::Helical180LI(int interpolation_view)
204 if (m_geometry == Scanner::GEOMETRY_INVALID)
206 std::cerr << "Invalid geometry " << m_geometry << std::endl;
209 else if (m_geometry == Scanner::GEOMETRY_PARALLEL)
211 std::cerr << "Helical 180LI not yet implemented for PARALLEL geometry"
215 else if (m_geometry == Scanner::GEOMETRY_EQUILINEAR)
217 std::cerr << "Helical 180LI not yet implemented for EQUILINEAR geometry"
221 else if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR)
223 return Helical180LI_Equiangular(interpolation_view);
227 std::cerr << "Invalid geometry in projection data file" << m_geometry
233 Projections::Helical180LI_Equiangular(int interpView)
235 double dbeta = m_rotInc;
236 double dgamma = m_detInc;
237 double fanAngle = m_dFanBeamAngle;
240 // is there enough data in the data set? Should have 2(Pi+fanAngle)
242 if ( m_nView < static_cast<int>((2*( PI + fanAngle ) ) / dbeta) -1 ){
243 std::cerr << "Data set does not include 360 +2*FanBeamAngle views"
248 if (interpView < 0) // use default position at PI+fanAngle
250 interpView = static_cast<int> ((PI+fanAngle)/dbeta);
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)
259 std::cerr << "There isn't PI+fanAngle of data on either side of the requested interpolation view" << std::endl;
262 offsetView = interpView - static_cast<int>((PI+fanAngle)/dbeta);
265 int last_interp_view = static_cast<int> ((PI+fanAngle)/dbeta);
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++)
279 int last_acq_view = 2*last_interp_view;
280 for ( int iView = 0 ; iView <= last_acq_view; iView++) {
281 double beta = iView * dbeta;
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 )
292 else // (beta > PI+fanAngle)
294 //newbeta = beta +2*gamma - 180;
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);
303 //std::cout << beta << " "<< gamma << " " << newbeta << " " << newgamma <<" " << iView-offsetView << " " << iDet << " " << newiView << " " << newiDet << std::endl;
304 //std::cout << iView-offsetView << " " << iDet << " " << newiView << " " << newiDet << std::endl;
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)
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)
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)
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)
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)
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)
349 ; // outside region of interest
355 m_projData = newdetarray;
356 m_nView = last_interp_view+1;
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.
371 Projections::HalfScanFeather(void)
373 double dbeta = m_rotInc;
374 double dgamma = m_detInc;
375 double fanAngle = m_dFanBeamAngle;
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;
382 if (m_geometry == Scanner::GEOMETRY_INVALID) {
383 std::cerr << "Invalid geometry " << m_geometry << std::endl;
387 if (m_geometry == Scanner::GEOMETRY_PARALLEL) {
388 std::cerr << "HalfScanFeather not yet implemented for PARALLEL geometry"<< std::endl;
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
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);
404 DetectorValue* detval2 = (m_projData[iView2])->detValues();
405 DetectorValue* detval1 = (m_projData[iView1])->detValues();
407 detval1[iDet1] = detval2[iDet2] ;
409 double x, w1,w2,beta1, gamma1;
412 if ( beta1 <= (fanAngle - 2*gamma1) )
413 x = beta1 / ( fanAngle - 2*gamma1);
414 else if ( (fanAngle - 2*gamma1 <= beta1 ) && beta1 <= PI - 2*gamma1)
416 else if ( (PI - 2*gamma1 <= beta1 ) && ( beta1 <=PI + fanAngle) )
417 x = (PI +fanAngle - beta1)/(fanAngle + 2*gamma1);
419 std::cerr << "Shouldn't be here!"<< std::endl;
422 w1 = (3*x - 2*x*x)*x;
424 detval1[iDet1] *= w1;
425 detval2[iDet2] *= w2;
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;
446 Projections::newProjData (void)
449 sys_error(ERR_WARNING, "m_projData != NULL [newProjData]");
451 if (m_nView > 0 && m_nDet) {
452 m_projData = new DetectorArray* [m_nView];
454 for (int i = 0; i < m_nView; i++)
455 m_projData[i] = new DetectorArray (m_nDet);
461 * projections_free Free memory allocated to projections
464 * projections_free(proj)
465 * Projections& proj Projectionss to be deallocated
469 Projections::deleteProjData (void)
471 if (m_projData != NULL) {
472 for (int i = 0; i < m_nView; i++)
473 delete m_projData[i];
482 * Projections::headerWwrite Write data header for projections file
487 Projections::headerWrite (fnetorderstream& fs)
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;
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;
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);
539 m_headerSize = fs.tellp();
540 _hsize = m_headerSize;
542 fs.writeInt16 (_hsize);
550 * projections_read_header Read data header for projections file
554 Projections::headerRead (fnetorderstream& fs)
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;
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);
581 fs.readInt16 (_hour);
582 fs.readInt16 (_minute);
583 fs.readInt16 (_second);
584 fs.readInt16 (_remarksize);
587 sys_error (ERR_SEVERE, "Error reading header information , _remarksize=%d [projections_read_header]", _remarksize);
591 if (_signature != m_signature) {
592 sys_error (ERR_SEVERE, "File %s does not have a valid projection file signature", m_filename.c_str());
596 char* pszRemarkStorage = new char [_remarksize+1];
597 fs.read (pszRemarkStorage, _remarksize);
599 sys_error (ERR_SEVERE, "Error reading remark, _remarksize = %d", _remarksize);
602 pszRemarkStorage[_remarksize] = 0;
603 m_remark = pszRemarkStorage;
604 delete pszRemarkStorage;
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);
612 m_headerSize = _hsize;
616 m_calcTime = _calcTime;
617 m_rotStart = _rotStart;
619 m_detStart = _detStart;
621 m_dFocalLength = _focalLength;
622 m_dSourceDetectorLength = _sourceDetectorLength;
623 m_dViewDiameter = _viewDiameter;
624 m_dFanBeamAngle = _fanBeamAngle;
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);
641 Projections::read (const std::string& filename)
643 return read (filename.c_str());
648 Projections::read (const char* filename)
650 m_filename = filename;
652 frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary);
654 frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary); // | std::ios::nocreate);
660 if (! headerRead (fileRead))
666 for (int i = 0; i < m_nView; i++) {
667 if (! detarrayRead (fileRead, *m_projData[i], i))
677 Projections::copyViewData (const std::string& filename, std::ostream& os, int startView, int endView)
679 return copyViewData (filename.c_str(), os, startView, endView);
683 Projections::copyViewData (const char* const filename, std::ostream& os, int startView, int endView)
685 frnetorderstream is (filename, std::ios::in | std::ios::binary);
686 kuint16 sizeHeader, signature;
687 kuint32 _nView, _nDet;
691 sys_error (ERR_SEVERE, "Unable to read projection file %s", filename);
695 is.readInt16 (sizeHeader);
696 is.readInt16 (signature);
697 is.readInt32 (_nView);
698 is.readInt32 (_nDet);
702 if (signature != m_signature) {
703 sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
709 if (startView > nView - 1)
711 if (endView < 0 || endView > nView - 1)
714 if (startView > endView) { // swap if start > end
715 int tempView = endView;
717 startView = tempView;
720 int sizeView = 8 /* view_angle */ + 4 /* nDet */ + (4 * nDet);
721 unsigned char* pViewData = new unsigned char [sizeView];
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())
733 sys_error (ERR_SEVERE, "Error reading projection file");
735 sys_error (ERR_SEVERE, "Error writing projection file");
737 return (! (is.fail() | os.fail()));
741 Projections::copyHeader (const std::string& filename, std::ostream& os)
743 return copyHeader (filename.c_str(), os);
747 Projections::copyHeader (const char* const filename, std::ostream& os)
749 frnetorderstream is (filename, std::ios::in | std::ios::binary);
750 kuint16 sizeHeader, signature;
751 is.readInt16 (sizeHeader);
752 is.readInt16 (signature);
754 if (signature != m_signature) {
755 sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
759 unsigned char* pHdrData = new unsigned char [sizeHeader];
760 is.read (reinterpret_cast<char*>(pHdrData), sizeHeader);
762 sys_error (ERR_SEVERE, "Error reading header");
766 os.write (reinterpret_cast<char*>(pHdrData), sizeHeader);
768 sys_error (ERR_SEVERE, "Error writing header");
776 Projections::write (const std::string& filename)
778 return write (filename.c_str());
782 Projections::write (const char* filename)
784 frnetorderstream fs (filename, std::ios::out | std::ios::binary | std::ios::trunc | std::ios::ate);
785 m_filename = filename;
787 sys_error (ERR_SEVERE, "Error opening file %s for output [projections_create]", filename);
791 if (! headerWrite (fs))
794 if (m_projData != NULL) {
795 for (int i = 0; i < m_nView; i++) {
796 if (! detarrayWrite (fs, *m_projData[i], i))
809 * detarrayRead Read a Detector Array structure from the disk
812 * detarrayRead (proj, darray, view_num)
813 * DETARRAY *darray Detector array storage location to be filled
814 * int view_num View number to read
818 Projections::detarrayRead (fnetorderstream& fs, DetectorArray& darray, const int iview)
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();
828 fs.seekg (start_data);
830 fs.readFloat64 (view_angle);
832 darray.setViewAngle (view_angle);
833 // darray.setNDet ( nDet);
835 for (unsigned int i = 0; i < nDet; i++) {
837 fs.readFloat32 (detval);
838 detval_ptr[i] = detval;
848 * detarrayWrite Write detector array data to the disk
851 * detarrayWrite (darray, view_num)
852 * DETARRAY *darray Detector array data to be written
853 * int view_num View number to write
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
862 Projections::detarrayWrite (fnetorderstream& fs, const DetectorArray& darray, const int iview)
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();
872 fs.seekp (start_data);
874 sys_error (ERR_SEVERE, "Error seeking detectory array [detarrayWrite]");
878 fs.writeFloat64 (view_angle);
879 fs.writeInt32 (nDet);
881 for (unsigned int i = 0; i < nDet; i++) {
882 kfloat32 detval = detval_ptr[i];
883 fs.writeFloat32 (detval);
893 * printProjectionData Print projections data
896 * printProjectionData ()
900 Projections::printProjectionData ()
902 printProjectionData (0, nView() - 1);
906 Projections::printProjectionData (int startView, int endView)
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) {
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]);
936 Projections::printScanInfo (std::ostringstream& os) const
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";
954 Projections::convertPolar (ImageFile& rIF, int iInterpolationID)
956 unsigned int nx = rIF.nx();
957 unsigned int ny = rIF.ny();
958 ImageFileArray v = rIF.getArray();
959 ImageFileArray vImag = rIF.getImaginaryArray();
961 if (! v || nx == 0 || ny == 0)
964 Projections* pProj = this;
965 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
966 pProj = interpolateToParallel();
968 Array2d<double> adView (nx, ny);
969 Array2d<double> adDet (nx, ny);
970 double** ppdView = adView.getArray();
971 double** ppdDet = adDet.getArray();
973 std::complex<double>** ppcDetValue = new std::complex<double>* [pProj->m_nView];
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);
982 pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet, pProj->m_nDet, 1., pProj->m_detInc);
984 pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet,
985 pProj->m_nDet, iInterpolationID);
987 for (iView = 0; iView < pProj->m_nView; iView++)
988 delete [] ppcDetValue[iView];
989 delete [] ppcDetValue;
991 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
999 Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad)
1002 rIF.arrayDataClear();
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();
1012 if (! v || nx == 0 || ny == 0)
1015 Projections* pProj = this;
1016 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
1017 pProj = interpolateToParallel();
1019 int iInterpDet = nx;
1020 // int iInterpDet = pProj->m_nDet;
1021 int iNumInterpDetWithZeros = ProcessSignal::addZeropadFactor (iInterpDet, iZeropad);
1023 double dZeropadRatio = static_cast<double>(iNumInterpDetWithZeros) / static_cast<double>(iInterpDet);
1025 fftw_plan plan = fftw_create_plan (iNumInterpDetWithZeros, FFTW_FORWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
1027 fftw_complex* pcIn = new fftw_complex [iNumInterpDetWithZeros];
1028 std::complex<double>** ppcDetValue = new std::complex<double>* [pProj->m_nView];
1029 double dInterpScale = (pProj->m_nDet-1) / static_cast<double>(iInterpDet-1) / SQRT2;
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;
1036 // For each view, interpolate to nx length, shift to center at origin, and FFt transform
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) * dInterpScale;
1046 Fourier::shuffleFourierToNaturalOrder (pcIn, iInterpDet);
1047 if (iZerosAdded > 0) {
1048 for (int iDet1 = iMidPoint; iDet1 < iInterpDet; iDet1++)
1049 pcIn[iDet1+iZerosAdded] = pcIn[iDet1];
1050 for (int iDet2 = iMidPoint; iDet2 < iMidPoint + iZerosAdded; iDet2++)
1051 pcIn[iDet2].re = pcIn[iDet2].im = 0;
1054 fftw_one (plan, pcIn, NULL);
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);
1061 Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], iNumInterpDetWithZeros);
1065 fftw_destroy_plan (plan);
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);
1074 pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet,
1075 iNumInterpDetWithZeros, iInterpolationID);
1077 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
1080 for (int i = 0; i < m_nView; i++)
1081 delete [] ppcDetValue[i];
1082 delete [] ppcDetValue;
1090 Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double** ppdView, double** ppdDet,
1091 int iNumDetWithZeros, double dZeropadRatio, double dDetInc)
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;
1102 xMin = (xMin - xCent) * dZeropadRatio + xCent;
1103 xMax = (xMax - xCent) * dZeropadRatio + xCent;
1104 yMin = (yMin - yCent) * dZeropadRatio + yCent;
1105 yMax = (yMax - yCent) * dZeropadRatio + yCent;
1107 double xInc = (xMax - xMin) / nx; // size of cells
1108 double yInc = (yMax - yMin) / ny;
1110 // +1 is correct for frequency data, ndet-1 is correct for projections
1111 int iDetCenter = (iNumDetWithZeros - 1) / 2; // index refering to L=0 projection
1112 if (isEven (iNumDetWithZeros))
1113 iDetCenter = (iNumDetWithZeros + 1) / 2;
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);
1130 ppdView[ix][iy] = (phi - m_rotStart) / m_rotInc;
1131 ppdDet[ix][iy] = (r / dDetInc) + iDetCenter;
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)
1141 typedef std::complex<double> complexValue;
1143 BilinearInterpolator<complexValue>* pBilinear = NULL;
1144 if (iInterpolationID == POLAR_INTERP_BILINEAR)
1145 pBilinear = new BilinearInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
1147 BicubicPolyInterpolator<complexValue>* pBicubic;
1148 if (iInterpolationID == POLAR_INTERP_BICUBIC)
1149 pBicubic = new BicubicPolyInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
1151 for (unsigned int ix = 0; ix < ny; ix++) {
1152 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]);
1157 if (iView == nView) {
1159 iDet = m_nDet - iDet;
1161 if (iDet >= 0 && iDet < nDetWithZeros && iView >= 0 && iView < nView) {
1162 v[ix][iy] = ppcDetValue[iView][iDet].real();
1164 vImag[ix][iy] = ppcDetValue[iView][iDet].imag();
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();
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();
1177 vImag[ix][iy] = vInterp.imag();
1184 Projections::initFromSomatomAR_STAR (int iNViews, int iNDets, unsigned char* pData, unsigned long lDataLength)
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);
1195 m_dViewDiameter = sin (m_dFanBeamAngle / 2) * m_dFocalLength * 2;
1197 if (! ((iNViews == 750 && lDataLength == 1560000L) || (iNViews == 950 && lDataLength == 1976000L)
1198 || (iNViews == 1500 && lDataLength == 3120000)))
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);
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);
1212 p = pArgBase+20; SwapBytes4IfLittleEndian (p);
1213 long lEscale = *reinterpret_cast<long*>(p);
1215 p = pArgBase+28; SwapBytes4IfLittleEndian (p);
1216 // long lTime = *reinterpret_cast<long*>(p);
1218 p = pArgBase + 4; SwapBytes4IfLittleEndian (p);
1219 double dAlpha = *reinterpret_cast<float*>(p) + HALFPI;
1221 p = pArgBase+12; SwapBytes4IfLittleEndian (p);
1222 // double dAlign = *reinterpret_cast<float*>(p);
1224 p = pArgBase + 16; SwapBytes4IfLittleEndian (p);
1225 // double dMaxValue = *reinterpret_cast<float*>(p);
1227 DetectorArray& detArray = getDetectorArray (iv);
1228 detArray.setViewAngle (dAlpha);
1229 DetectorValue* detval = detArray.detValues();
1231 double dViewScale = 1. / (2294.4871 * ::pow (2.0, -lEscale));
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
1237 detval[id] = iV * dViewScale * pdCosScale[id];
1241 for (int k = iNDets - 2; k >= 0; k--)
1242 detval[k+1] = detval[k];
1252 Projections::interpolateToParallel () const
1254 if (m_geometry == Scanner::GEOMETRY_PARALLEL)
1255 return const_cast<Projections*>(this);
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);
1273 pProjNew->m_rotStart = 0;
1274 #ifdef CONVERT_PARALLEL_PI
1275 pProjNew->m_rotInc = PI / nView;;
1277 pProjNew->m_rotInc = TWOPI / nView;
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);
1284 ParallelRaysums parallel (this, ParallelRaysums::THETA_RANGE_NORMALIZE_TO_TWOPI);
1286 double* pdThetaValuesForT = new double [pProjNew->nView()];
1287 double* pdRaysumsForT = new double [pProjNew->nView()];
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);
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);
1302 delete pdThetaValuesForT;
1303 delete pdRaysumsForT;
1305 // interpolate to evenly space t (detectors)
1306 double* pdOriginalDetPositions = new double [pProjNew->nDet()];
1307 parallel.getDetPositions (pdOriginalDetPositions);
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);
1316 for (int i = 0; i < pProjNew->nDet(); i++)
1317 pdDetValueCopy[i] = detValues[i];
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);
1325 delete pdDetValueCopy;
1326 delete pdOriginalDetPositions;
1332 ///////////////////////////////////////////////////////////////////////////////
1334 // Class ParallelRaysums
1336 // Used for converting divergent beam raysums into Parallel raysums
1338 ///////////////////////////////////////////////////////////////////////////////
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)
1344 int iGeometry = pProjections->geometry();
1345 double dDetInc = pProjections->detInc();
1346 double dDetStart = pProjections->detStart();
1347 double dFocalLength = pProjections->focalLength();
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;
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();
1360 double dDetPos = dDetStart;
1361 for (int iD = 0; iD < m_iNumDet; iD++) {
1362 ParallelRaysumCoordinate* pC = m_vecpCoordinates[iCoordinate++];
1364 if (iGeometry == Scanner::GEOMETRY_PARALLEL) {
1365 pC->m_dTheta = dViewAngle;
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);
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);
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) {
1381 pC->m_dT = -pC->m_dT;
1384 pC->m_dRaysum = detValues[iD];
1390 ParallelRaysums::~ParallelRaysums()
1392 delete m_pCoordinates;
1395 ParallelRaysums::CoordinateContainer&
1396 ParallelRaysums::getSortedByTheta()
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);
1405 return m_vecpSortedByTheta;
1408 ParallelRaysums::CoordinateContainer&
1409 ParallelRaysums::getSortedByT()
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);
1418 return m_vecpSortedByT;
1423 ParallelRaysums::getLimits (double* dMinT, double* dMaxT, double* dMinTheta, double* dMaxTheta) const
1425 if (m_iNumCoordinates <= 0)
1428 *dMinT = *dMaxT = m_vecpCoordinates[0]->m_dT;
1429 *dMinTheta = *dMaxTheta = m_vecpCoordinates[0]->m_dTheta;
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;
1437 else if (dT > *dMaxT)
1440 if (dTheta < *dMinTheta)
1441 *dMinTheta = dTheta;
1442 else if (dTheta > *dMaxTheta)
1443 *dMaxTheta = dTheta;
1448 ParallelRaysums::getThetaAndRaysumsForT (int iTheta, double* pTheta, double* pRaysum)
1450 const CoordinateContainer& coordsT = getSortedByT();
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;
1461 ParallelRaysums::getDetPositions (double* pdDetPos)
1463 const CoordinateContainer& coordsT = getSortedByT();
1466 for (int i = 0; i < m_iNumDet; i++) {
1467 pdDetPos[i] = coordsT[iPos]->m_dT;