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.77 2002/05/28 18:43:16 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 ******************************************************************************/
29 const kuint16 Projections::m_signature = ('P'*256 + 'J');
31 const int Projections::POLAR_INTERP_INVALID = -1;
32 const int Projections::POLAR_INTERP_NEAREST = 0;
33 const int Projections::POLAR_INTERP_BILINEAR = 1;
34 const int Projections::POLAR_INTERP_BICUBIC = 2;
36 const char* const Projections::s_aszInterpName[] =
43 const char* const Projections::s_aszInterpTitle[] =
50 const int Projections::s_iInterpCount = sizeof(s_aszInterpName) / sizeof(char*);
55 * Projections Constructor for projections matrix storage
58 * proj = projections_create (filename, nView, nDet)
59 * Projections& proj Allocated projections structure & matrix
60 * int nView Number of rotated view
61 * int nDet Number of detectors
65 Projections::Projections (const Scanner& scanner)
68 initFromScanner (scanner);
72 Projections::Projections (const int nView, const int nDet)
78 Projections::Projections (void)
84 Projections::~Projections (void)
90 Projections::convertInterpNameToID (const char* const interpName)
92 int interpID = POLAR_INTERP_INVALID;
94 for (int i = 0; i < s_iInterpCount; i++)
95 if (strcasecmp (interpName, s_aszInterpName[i]) == 0) {
104 Projections::convertInterpIDToName (const int interpID)
106 static const char *interpName = "";
108 if (interpID >= 0 && interpID < s_iInterpCount)
109 return (s_aszInterpName[interpID]);
115 Projections::convertInterpIDToTitle (const int interpID)
117 static const char *interpTitle = "";
119 if (interpID >= 0 && interpID < s_iInterpCount)
120 return (s_aszInterpTitle[interpID]);
122 return (interpTitle);
128 Projections::init (const int nView, const int nDet)
130 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
135 time_t t = time (NULL);
136 tm* lt = localtime (&t);
137 m_year = lt->tm_year;
138 m_month = lt->tm_mon;
140 m_hour = lt->tm_hour;
141 m_minute = lt->tm_min;
142 m_second = lt->tm_sec;
146 Projections::initFromScanner (const Scanner& scanner)
148 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
150 init (scanner.nView(), scanner.nDet());
152 m_rotInc = scanner.rotInc();
153 m_detInc = scanner.detInc();
154 m_detStart = scanner.detStart();
155 m_geometry = scanner.geometry();
156 m_dFocalLength = scanner.focalLength();
157 m_dSourceDetectorLength = scanner.sourceDetectorLength();
158 m_dViewDiameter = scanner.viewDiameter();
159 m_rotStart = scanner.offsetView()*scanner.rotInc();
160 m_dFanBeamAngle = scanner.fanBeamAngle();
164 Projections::setNView (int nView) // used by MPI to reduce # of views
167 init (nView, m_nDet);
170 // Helical 180 Linear Interpolation.
171 // This member function takes a set of helical scan projections and
172 // performs a linear interpolation between pairs of complementary rays
173 // to produce a single projection data set approximating what would be
174 // measured at a single axial plane.
175 // Complementary rays are rays which traverse the same path through the
176 // phantom in opposite directions.
178 // For parallel beam geometry, a ray with a given gantry angle beta and a
179 // detector iDet will have a complementary ray at beta + pi and nDet-iDet
181 // For equiangular or equilinear beam geometry the complementary ray to
182 // gantry angle beta and fan-beam angle gamma is at
183 // beta-hat = beta +2*gamma + pi, and gamma-hat = -gamma.
184 // Note that beta-hat - beta depends on gamma and is not constant.
186 // The algorithm used here is from Crawford and King, Med. Phys. 17(6)
187 // 1990 p967; what they called method "C", CSH-HH. It uses interpolation only
188 // between pairs of complementary rays on either side of an image plane.
189 // Input data must sample gantry angles from zero to
190 // (2*pi + 2* fan-beam-angle). The data set produced contains gantry
191 // angles from 0 to Pi+fan-beam-angle. This is a "halfscan" data set,
192 // which still contains redundant data, and can be used with a half scan
193 // reconstruction to produce an image.
194 // In this particular implementation a lower triangle from (beta,gamma) =
195 // (0,-fanAngle/2)->(2*fanAngle,-fanAngle/2)->(0,fanAngle/2) contains
196 // zeros, but is actually redundant with data contained in the region
197 // (pi+fanAngle,-fanAngle/2)->(pi+fanAngle, fanAngle/2) ->(pi-fanAngle,
201 Projections::Helical180LI(int interpolation_view)
203 if (m_geometry == Scanner::GEOMETRY_INVALID)
205 std::cerr << "Invalid geometry " << m_geometry << std::endl;
208 else if (m_geometry == Scanner::GEOMETRY_PARALLEL)
210 std::cerr << "Helical 180LI not yet implemented for PARALLEL geometry"
214 else if (m_geometry == Scanner::GEOMETRY_EQUILINEAR)
216 std::cerr << "Helical 180LI not yet implemented for EQUILINEAR geometry"
220 else if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR)
222 return Helical180LI_Equiangular(interpolation_view);
226 std::cerr << "Invalid geometry in projection data file" << m_geometry
232 Projections::Helical180LI_Equiangular(int interpView)
234 double dbeta = m_rotInc;
235 double dgamma = m_detInc;
236 double fanAngle = m_dFanBeamAngle;
239 // is there enough data in the data set? Should have 2(Pi+fanAngle)
241 if ( m_nView < static_cast<int>((2*( PI + fanAngle ) ) / dbeta) -1 ){
242 std::cerr << "Data set does not include 360 +2*FanBeamAngle views"
247 if (interpView < 0) // use default position at PI+fanAngle
249 interpView = static_cast<int> ((PI+fanAngle)/dbeta);
253 // check if there is PI+fanAngle data on either side of the
254 // of the specified image plane
255 if ( interpView*dbeta < PI+fanAngle ||
256 interpView*dbeta + PI + fanAngle > m_nView*dbeta)
258 std::cerr << "There isn't PI+fanAngle of data on either side of the requested interpolation view" << std::endl;
261 offsetView = interpView - static_cast<int>((PI+fanAngle)/dbeta);
264 int last_interp_view = static_cast<int> ((PI+fanAngle)/dbeta);
267 // make a new array for data...
268 class DetectorArray ** newdetarray = new DetectorArray * [last_interp_view+1];
269 for ( int i=0 ; i <= last_interp_view ; i++ ){
270 newdetarray[i] = new DetectorArray (m_nDet);
271 newdetarray[i]->setViewAngle((i+offsetView)*dbeta);
272 DetectorValue* newdetval = (newdetarray[i])->detValues();
273 // and initialize the data to zero
274 for (int j=0; j < m_nDet; j++)
278 int last_acq_view = 2*last_interp_view;
279 for ( int iView = 0 ; iView <= last_acq_view; iView++) {
280 double beta = iView * dbeta;
282 for ( int iDet = 0; iDet < m_nDet; iDet++) {
283 double gamma = (iDet -(m_nDet-1)/2)* dgamma ;
284 int newiView, newiDet;
285 if (beta < PI+fanAngle) { //if (PI +fanAngle - beta > dbeta )
291 else // (beta > PI+fanAngle)
293 //newbeta = beta +2*gamma - 180;
295 newiDet = -iDet + (m_nDet -1);
296 // newiView = nearest<int>((beta + 2*gamma - PI)/dbeta);
297 //newiView = static_cast<int>(( (iView*dbeta) + 2*(iDet-(m_nDet-1)/2)*dgamma - PI)/dbeta);
298 newiView = nearest<int>(( (iView*dbeta) + 2*(iDet-(m_nDet-1)/2)*dgamma - PI)/dbeta);
302 //std::cout << beta << " "<< gamma << " " << newbeta << " " << newgamma <<" " << iView-offsetView << " " << iDet << " " << newiView << " " << newiDet << std::endl;
303 //std::cout << iView-offsetView << " " << iDet << " " << newiView << " " << newiDet << std::endl;
306 if ( ( beta > fanAngle - 2*gamma)
307 && ( beta < 2*PI + fanAngle -2*gamma) )
308 { // not in region 1 or 8
309 DetectorValue* detval = (m_projData[iView+offsetView])->detValues();
310 DetectorValue* newdetval = (newdetarray[newiView])->detValues();
311 if ( beta > fanAngle - 2*gamma
312 && beta <= 2*fanAngle ) { // in region 2
313 newdetval[newiDet] +=
314 (beta +2*gamma - fanAngle)/(PI+2*gamma)
316 } else if ( beta > 2*fanAngle
317 && beta <= PI - 2*gamma) { // in region 3
318 newdetval[newiDet] +=
319 (beta +2*gamma - fanAngle)/(PI+2*gamma)
322 else if ( beta > PI -2*gamma
323 && beta <= PI + fanAngle ) { // in region 4
324 newdetval[newiDet] +=
325 (beta +2*gamma - fanAngle)/(PI+2*gamma)
328 else if ( beta > PI + fanAngle
329 && beta <= PI +2*fanAngle -2*gamma) { // in region 5
330 newdetval[newiDet] +=
331 (2*PI - beta - 2*gamma + fanAngle)/(PI-2*gamma)
334 else if ( beta > PI +2*fanAngle -2*gamma
335 && beta <= 2*PI) { // in region 6
336 newdetval[newiDet] +=
337 (2*PI - beta - 2*gamma + fanAngle)/(PI-2*gamma)
340 else if ( beta > 2*PI
341 && beta <= 2*PI + fanAngle -2*gamma){ // in region 7
342 newdetval[newiDet] +=
343 (2*PI - beta -2*gamma + fanAngle)/(PI-2*gamma)
348 ; // outside region of interest
354 m_projData = newdetarray;
355 m_nView = last_interp_view+1;
360 // A HalfScan Projection Data Set for equiangular geometry,
361 // covering gantry angles from 0 to pi+fanBeamAngle
362 // and fan angle gamma from -fanBeamAngle/2 to fanBeamAngle/2
363 // contains redundant information. If one copy of this data is left as
364 // zero, (as in the Helical180LI routine above) overweighting is avoided,
365 // but the discontinuity in the data introduces ringing in the image.
366 // This routine makes a copy of the data and applies a weighting to avoid
367 // over-representation, as given in Appendix C of Crawford and King, Med
368 // Phys 17 1990, p967.
370 Projections::HalfScanFeather(void)
372 double dbeta = m_rotInc;
373 double dgamma = m_detInc;
374 double fanAngle = m_dFanBeamAngle;
376 // is there enough data?
377 if ( m_nView != static_cast<int>(( PI+fanAngle ) / dbeta) +1 ){
378 std::cerr << "Data set does seem have enough data to be a halfscan data set" << std::endl;
381 if (m_geometry == Scanner::GEOMETRY_INVALID) {
382 std::cerr << "Invalid geometry " << m_geometry << std::endl;
386 if (m_geometry == Scanner::GEOMETRY_PARALLEL) {
387 std::cerr << "HalfScanFeather not yet implemented for PARALLEL geometry"<< std::endl;
391 for ( int iView2 = 0 ; iView2 < m_nView; iView2++) {
392 double beta2 = iView2 * dbeta;
393 for ( int iDet2 = 0; iDet2 < m_nDet; iDet2++) {
394 double gamma2 = (iDet2 -(m_nDet-1)/2)* dgamma ;
395 if ( ( beta2 >= PI - 2*gamma2) ) { // in redundant data region
397 iDet1 = (m_nDet -1) - iDet2;
398 //iView1 = nearest<int>((beta2 + 2*gamma2 - PI)/dbeta);
399 iView1 = nearest<int>(( (iView2*dbeta)
400 + 2*(iDet2-(m_nDet-1)/2)*dgamma - PI)/dbeta);
403 DetectorValue* detval2 = (m_projData[iView2])->detValues();
404 DetectorValue* detval1 = (m_projData[iView1])->detValues();
406 detval1[iDet1] = detval2[iDet2] ;
408 double x, w1,w2,beta1, gamma1;
411 if ( beta1 <= (fanAngle - 2*gamma1) )
412 x = beta1 / ( fanAngle - 2*gamma1);
413 else if ( (fanAngle - 2*gamma1 <= beta1 ) && beta1 <= PI - 2*gamma1)
415 else if ( (PI - 2*gamma1 <= beta1 ) && ( beta1 <=PI + fanAngle) )
416 x = (PI +fanAngle - beta1)/(fanAngle + 2*gamma1);
418 std::cerr << "Shouldn't be here!"<< std::endl;
421 w1 = (3*x - 2*x*x)*x;
423 detval1[iDet1] *= w1;
424 detval2[iDet2] *= w2;
429 // heuristic scaling, why this factor?
430 double scalefactor = m_nView * m_rotInc / PI;
431 for ( int iView = 0 ; iView < m_nView; iView++) {
432 DetectorValue* detval = (m_projData[iView])->detValues();
433 for ( int iDet = 0; iDet < m_nDet; iDet++) {
434 detval[iDet] *= scalefactor;
445 Projections::newProjData (void)
448 sys_error(ERR_WARNING, "m_projData != NULL [newProjData]");
450 if (m_nView > 0 && m_nDet) {
451 m_projData = new DetectorArray* [m_nView];
453 for (int i = 0; i < m_nView; i++)
454 m_projData[i] = new DetectorArray (m_nDet);
460 * projections_free Free memory allocated to projections
463 * projections_free(proj)
464 * Projections& proj Projectionss to be deallocated
468 Projections::deleteProjData (void)
470 if (m_projData != NULL) {
471 for (int i = 0; i < m_nView; i++)
472 delete m_projData[i];
481 * Projections::headerWwrite Write data header for projections file
486 Projections::headerWrite (fnetorderstream& fs)
488 kuint16 _hsize = m_headerSize;
489 kuint16 _signature = m_signature;
490 kuint32 _nView = m_nView;
491 kuint32 _nDet = m_nDet;
492 kuint32 _geom = m_geometry;
493 kuint16 _remarksize = m_remark.length();
494 kuint16 _year = m_year;
495 kuint16 _month = m_month;
496 kuint16 _day = m_day;
497 kuint16 _hour = m_hour;
498 kuint16 _minute = m_minute;
499 kuint16 _second = m_second;
501 kfloat64 _calcTime = m_calcTime;
502 kfloat64 _rotStart = m_rotStart;
503 kfloat64 _rotInc = m_rotInc;
504 kfloat64 _detStart = m_detStart;
505 kfloat64 _detInc = m_detInc;
506 kfloat64 _viewDiameter = m_dViewDiameter;
507 kfloat64 _focalLength = m_dFocalLength;
508 kfloat64 _sourceDetectorLength = m_dSourceDetectorLength;
509 kfloat64 _fanBeamAngle = m_dFanBeamAngle;
515 fs.writeInt16 (_hsize);
516 fs.writeInt16 (_signature);
517 fs.writeInt32 (_nView);
518 fs.writeInt32 (_nDet);
519 fs.writeInt32 (_geom);
520 fs.writeFloat64 (_calcTime);
521 fs.writeFloat64 (_rotStart);
522 fs.writeFloat64 (_rotInc);
523 fs.writeFloat64 (_detStart);
524 fs.writeFloat64 (_detInc);
525 fs.writeFloat64 (_viewDiameter);
526 fs.writeFloat64 (_focalLength);
527 fs.writeFloat64 (_sourceDetectorLength);
528 fs.writeFloat64 (_fanBeamAngle);
529 fs.writeInt16 (_year);
530 fs.writeInt16 (_month);
531 fs.writeInt16 (_day);
532 fs.writeInt16 (_hour);
533 fs.writeInt16 (_minute);
534 fs.writeInt16 (_second);
535 fs.writeInt16 (_remarksize);
536 fs.write (m_remark.c_str(), _remarksize);
538 m_headerSize = fs.tellp();
539 _hsize = m_headerSize;
541 fs.writeInt16 (_hsize);
549 * projections_read_header Read data header for projections file
553 Projections::headerRead (fnetorderstream& fs)
555 kuint16 _hsize, _signature, _year, _month, _day, _hour, _minute, _second, _remarksize = 0;
556 kuint32 _nView, _nDet, _geom;
557 kfloat64 _calcTime, _rotStart, _rotInc, _detStart, _detInc, _focalLength, _sourceDetectorLength, _viewDiameter, _fanBeamAngle;
563 fs.readInt16 (_hsize);
564 fs.readInt16 (_signature);
565 fs.readInt32 (_nView);
566 fs.readInt32 (_nDet);
567 fs.readInt32 (_geom);
568 fs.readFloat64 (_calcTime);
569 fs.readFloat64 (_rotStart);
570 fs.readFloat64 (_rotInc);
571 fs.readFloat64 (_detStart);
572 fs.readFloat64 (_detInc);
573 fs.readFloat64 (_viewDiameter);
574 fs.readFloat64 (_focalLength);
575 fs.readFloat64 (_sourceDetectorLength);
576 fs.readFloat64 (_fanBeamAngle);
577 fs.readInt16 (_year);
578 fs.readInt16 (_month);
580 fs.readInt16 (_hour);
581 fs.readInt16 (_minute);
582 fs.readInt16 (_second);
583 fs.readInt16 (_remarksize);
586 sys_error (ERR_SEVERE, "Error reading header information , _remarksize=%d [projections_read_header]", _remarksize);
590 if (_signature != m_signature) {
591 sys_error (ERR_SEVERE, "File %s does not have a valid projection file signature", m_filename.c_str());
595 char* pszRemarkStorage = new char [_remarksize+1];
596 fs.read (pszRemarkStorage, _remarksize);
598 sys_error (ERR_SEVERE, "Error reading remark, _remarksize = %d", _remarksize);
601 pszRemarkStorage[_remarksize] = 0;
602 m_remark = pszRemarkStorage;
603 delete pszRemarkStorage;
605 off_t _hsizeread = fs.tellg();
606 if (!fs || _hsizeread != _hsize) {
607 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);
611 m_headerSize = _hsize;
615 m_calcTime = _calcTime;
616 m_rotStart = _rotStart;
618 m_detStart = _detStart;
620 m_dFocalLength = _focalLength;
621 m_dSourceDetectorLength = _sourceDetectorLength;
622 m_dViewDiameter = _viewDiameter;
623 m_dFanBeamAngle = _fanBeamAngle;
631 m_label.setLabelType (Array2dFileLabel::L_HISTORY);
632 m_label.setLabelString (m_remark);
633 m_label.setCalcTime (m_calcTime);
634 m_label.setDateTime (m_year, m_month, m_day, m_hour, m_minute, m_second);
640 Projections::read (const std::string& filename)
642 return read (filename.c_str());
647 Projections::read (const char* filename)
649 m_filename = filename;
651 frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary);
653 frnetorderstream fileRead (m_filename.c_str(), std::ios::in | std::ios::binary); // | std::ios::nocreate);
659 if (! headerRead (fileRead))
665 for (int i = 0; i < m_nView; i++) {
666 if (! detarrayRead (fileRead, *m_projData[i], i))
676 Projections::copyViewData (const std::string& filename, std::ostream& os, int startView, int endView)
678 return copyViewData (filename.c_str(), os, startView, endView);
682 Projections::copyViewData (const char* const filename, std::ostream& os, int startView, int endView)
684 frnetorderstream is (filename, std::ios::in | std::ios::binary);
685 kuint16 sizeHeader, signature;
686 kuint32 _nView, _nDet;
690 sys_error (ERR_SEVERE, "Unable to read projection file %s", filename);
694 is.readInt16 (sizeHeader);
695 is.readInt16 (signature);
696 is.readInt32 (_nView);
697 is.readInt32 (_nDet);
701 if (signature != m_signature) {
702 sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
708 if (startView > nView - 1)
710 if (endView < 0 || endView > nView - 1)
713 if (startView > endView) { // swap if start > end
714 int tempView = endView;
716 startView = tempView;
719 int sizeView = 8 /* view_angle */ + 4 /* nDet */ + (4 * nDet);
720 unsigned char* pViewData = new unsigned char [sizeView];
722 for (int i = startView; i <= endView; i++) {
723 is.seekg (sizeHeader + i * sizeView);
724 is.read (reinterpret_cast<char*>(pViewData), sizeView);
725 os.write (reinterpret_cast<char*>(pViewData), sizeView);
726 if (is.fail() || os.fail())
732 sys_error (ERR_SEVERE, "Error reading projection file");
734 sys_error (ERR_SEVERE, "Error writing projection file");
736 return (! (is.fail() | os.fail()));
740 Projections::copyHeader (const std::string& filename, std::ostream& os)
742 return copyHeader (filename.c_str(), os);
746 Projections::copyHeader (const char* const filename, std::ostream& os)
748 frnetorderstream is (filename, std::ios::in | std::ios::binary);
749 kuint16 sizeHeader, signature;
750 is.readInt16 (sizeHeader);
751 is.readInt16 (signature);
753 if (signature != m_signature) {
754 sys_error (ERR_SEVERE, "Illegal signature in projection file %s", filename);
758 unsigned char* pHdrData = new unsigned char [sizeHeader];
759 is.read (reinterpret_cast<char*>(pHdrData), sizeHeader);
761 sys_error (ERR_SEVERE, "Error reading header");
765 os.write (reinterpret_cast<char*>(pHdrData), sizeHeader);
767 sys_error (ERR_SEVERE, "Error writing header");
775 Projections::write (const std::string& filename)
777 return write (filename.c_str());
781 Projections::write (const char* filename)
783 frnetorderstream fs (filename, std::ios::out | std::ios::binary | std::ios::trunc | std::ios::ate);
784 m_filename = filename;
786 sys_error (ERR_SEVERE, "Error opening file %s for output [projections_create]", filename);
790 if (! headerWrite (fs))
793 if (m_projData != NULL) {
794 for (int i = 0; i < m_nView; i++) {
795 if (! detarrayWrite (fs, *m_projData[i], i))
808 * detarrayRead Read a Detector Array structure from the disk
811 * detarrayRead (proj, darray, view_num)
812 * DETARRAY *darray Detector array storage location to be filled
813 * int view_num View number to read
817 Projections::detarrayRead (fnetorderstream& fs, DetectorArray& darray, const int iview)
819 const int detval_bytes = darray.nDet() * sizeof(kfloat32);
820 const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
821 const int view_bytes = detheader_bytes + detval_bytes;
822 const off_t start_data = m_headerSize + (iview * view_bytes);
823 DetectorValue* detval_ptr = darray.detValues();
827 fs.seekg (start_data);
829 fs.readFloat64 (view_angle);
831 darray.setViewAngle (view_angle);
832 // darray.setNDet ( nDet);
834 for (unsigned int i = 0; i < nDet; i++) {
836 fs.readFloat32 (detval);
837 detval_ptr[i] = detval;
847 * detarrayWrite Write detector array data to the disk
850 * detarrayWrite (darray, view_num)
851 * DETARRAY *darray Detector array data to be written
852 * int view_num View number to write
855 * This routine writes the detarray data from the disk sequentially to
856 * the file that was opened with open_projections(). Data is written in
861 Projections::detarrayWrite (fnetorderstream& fs, const DetectorArray& darray, const int iview)
863 const int detval_bytes = darray.nDet() * sizeof(float);
864 const int detheader_bytes = sizeof(kfloat64) /* view_angle */ + sizeof(kint32) /* nDet */;
865 const int view_bytes = detheader_bytes + detval_bytes;
866 const off_t start_data = m_headerSize + (iview * view_bytes);
867 const DetectorValue* const detval_ptr = darray.detValues();
868 kfloat64 view_angle = darray.viewAngle();
869 kuint32 nDet = darray.nDet();
871 fs.seekp (start_data);
873 sys_error (ERR_SEVERE, "Error seeking detectory array [detarrayWrite]");
877 fs.writeFloat64 (view_angle);
878 fs.writeInt32 (nDet);
880 for (unsigned int i = 0; i < nDet; i++) {
881 kfloat32 detval = detval_ptr[i];
882 fs.writeFloat32 (detval);
892 * printProjectionData Print projections data
895 * printProjectionData ()
899 Projections::printProjectionData ()
901 printProjectionData (0, nView() - 1);
905 Projections::printProjectionData (int startView, int endView)
907 printf("Projections Data\n\n");
908 printf("Description: %s\n", m_remark.c_str());
909 printf("Geometry: %s\n", Scanner::convertGeometryIDToName (m_geometry));
910 printf("nView = %8d nDet = %8d\n", m_nView, m_nDet);
911 printf("focalLength = %8.4f ViewDiameter = %8.4f\n", m_dFocalLength, m_dViewDiameter);
912 printf("fanBeamAngle= %8.4f SourceDetector = %8.4f\n", convertRadiansToDegrees(m_dFanBeamAngle), m_dSourceDetectorLength);
913 printf("rotStart = %8.4f rotInc = %8.4f\n", m_rotStart, m_rotInc);
914 printf("detStart = %8.4f detInc = %8.4f\n", m_detStart, m_detInc);
915 if (m_projData != NULL) {
919 endView = m_nView - 1;
920 if (startView > m_nView - 1)
921 startView = m_nView - 1;
922 if (endView > m_nView - 1)
923 endView = m_nView - 1;
924 for (int ir = startView; ir <= endView - 1; ir++) {
925 printf("View %d: angle %f\n", ir, m_projData[ir]->viewAngle());
926 DetectorValue* detval = m_projData[ir]->detValues();
927 for (int id = 0; id < m_projData[ir]->nDet(); id++)
928 printf("%8.4f ", detval[id]);
935 Projections::printScanInfo (std::ostringstream& os) const
937 os << "Number of detectors: " << m_nDet << "\n";
938 os << "Number of views: " << m_nView<< "\n";
939 os << "Description: " << m_remark.c_str()<< "\n";
940 os << "Geometry: " << Scanner::convertGeometryIDToName (m_geometry)<< "\n";
941 os << "Focal Length: " << m_dFocalLength<< "\n";
942 os << "Source Detector Length: " << m_dSourceDetectorLength << "\n";
943 os << "View Diameter: " << m_dViewDiameter<< "\n";
944 os << "Fan Beam Angle: " << convertRadiansToDegrees(m_dFanBeamAngle) << "\n";
945 os << "detStart: " << m_detStart<< "\n";
946 os << "detInc: " << m_detInc<< "\n";
947 os << "rotStart: " << m_rotStart<< "\n";
948 os << "rotInc: " << m_rotInc<< "\n";
953 Projections::convertPolar (ImageFile& rIF, int iInterpolationID)
955 unsigned int nx = rIF.nx();
956 unsigned int ny = rIF.ny();
957 ImageFileArray v = rIF.getArray();
958 ImageFileArray vImag = rIF.getImaginaryArray();
960 if (! v || nx == 0 || ny == 0)
963 Projections* pProj = this;
964 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
965 pProj = interpolateToParallel();
967 Array2d<double> adView (nx, ny);
968 Array2d<double> adDet (nx, ny);
969 double** ppdView = adView.getArray();
970 double** ppdDet = adDet.getArray();
972 std::complex<double>** ppcDetValue = new std::complex<double>* [pProj->m_nView];
974 for (iView = 0; iView < pProj->m_nView; iView++) {
975 ppcDetValue[iView] = new std::complex<double> [pProj->m_nDet];
976 DetectorValue* detval = pProj->getDetectorArray (iView).detValues();
977 for (unsigned int iDet = 0; iDet < pProj->m_nDet; iDet++)
978 ppcDetValue[iView][iDet] = std::complex<double>(detval[iDet], 0);
981 pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet, pProj->m_nDet, 1., pProj->m_detInc);
983 pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet,
984 pProj->m_nDet, iInterpolationID);
986 for (iView = 0; iView < pProj->m_nView; iView++)
987 delete [] ppcDetValue[iView];
988 delete [] ppcDetValue;
990 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
998 Projections::convertFFTPolar (ImageFile& rIF, int iInterpolationID, int iZeropad)
1001 rIF.arrayDataClear();
1004 unsigned int nx = rIF.nx();
1005 unsigned int ny = rIF.ny();
1006 ImageFileArray v = rIF.getArray();
1007 if (! rIF.isComplex())
1008 rIF.convertRealToComplex();
1009 ImageFileArray vImag = rIF.getImaginaryArray();
1011 if (! v || nx == 0 || ny == 0)
1014 Projections* pProj = this;
1015 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
1016 pProj = interpolateToParallel();
1018 int iInterpDet = nx;
1019 // int iInterpDet = pProj->m_nDet;
1020 int iNumInterpDetWithZeros = ProcessSignal::addZeropadFactor (iInterpDet, iZeropad);
1022 double dZeropadRatio = static_cast<double>(iNumInterpDetWithZeros) / static_cast<double>(iInterpDet);
1024 fftw_plan plan = fftw_create_plan (iNumInterpDetWithZeros, FFTW_FORWARD, FFTW_IN_PLACE | FFTW_ESTIMATE | FFTW_USE_WISDOM);
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) / SQRT2;
1030 double dFFTScale = 1. / static_cast<double>(iInterpDet * iInterpDet);
1031 int iMidPoint = iInterpDet / 2;
1032 double dMidPoint = static_cast<double>(iInterpDet) / 2.;
1033 int iZerosAdded = iNumInterpDetWithZeros - iInterpDet;
1035 // For each view, interpolate to nx length, shift to center at origin, and FFt transform
1036 for (unsigned int iView = 0; iView < m_nView; iView++) {
1037 DetectorValue* detval = pProj->getDetectorArray(iView).detValues();
1038 LinearInterpolator<DetectorValue> projInterp (detval, pProj->m_nDet);
1039 for (unsigned int iDet = 0; iDet < iInterpDet; iDet++) {
1040 double dInterpPos = (m_nDet / 2.) + (iDet - dMidPoint) * dInterpScale;
1041 pcIn[iDet].re = projInterp.interpolate (dInterpPos) * dInterpScale;
1045 Fourier::shuffleFourierToNaturalOrder (pcIn, iInterpDet);
1046 if (iZerosAdded > 0) {
1047 for (unsigned int iDet1 = iMidPoint; iDet1 < iInterpDet; iDet1++)
1048 pcIn[iDet1+iZerosAdded] = pcIn[iDet1];
1049 for (unsigned int iDet2 = iMidPoint; iDet2 < iMidPoint + iZerosAdded; iDet2++)
1050 pcIn[iDet2].re = pcIn[iDet2].im = 0;
1053 fftw_one (plan, pcIn, NULL);
1055 ppcDetValue[iView] = new std::complex<double> [iNumInterpDetWithZeros];
1056 for (unsigned int iD = 0; iD < iNumInterpDetWithZeros; iD++) {
1057 ppcDetValue[iView][iD] = std::complex<double> (pcIn[iD].re * dFFTScale, pcIn[iD].im * dFFTScale);
1060 Fourier::shuffleFourierToNaturalOrder (ppcDetValue[iView], iNumInterpDetWithZeros);
1064 fftw_destroy_plan (plan);
1066 Array2d<double> adView (nx, ny);
1067 Array2d<double> adDet (nx, ny);
1068 double** ppdView = adView.getArray();
1069 double** ppdDet = adDet.getArray();
1070 pProj->calcArrayPolarCoordinates (nx, ny, ppdView, ppdDet, iNumInterpDetWithZeros, dZeropadRatio,
1071 pProj->m_detInc * dInterpScale);
1073 pProj->interpolatePolar (v, vImag, nx, ny, ppcDetValue, ppdView, ppdDet, pProj->m_nView, pProj->m_nDet,
1074 iNumInterpDetWithZeros, iInterpolationID);
1076 if (m_geometry == Scanner::GEOMETRY_EQUIANGULAR || m_geometry == Scanner::GEOMETRY_EQUILINEAR)
1079 for (int i = 0; i < m_nView; i++)
1080 delete [] ppcDetValue[i];
1081 delete [] ppcDetValue;
1089 Projections::calcArrayPolarCoordinates (unsigned int nx, unsigned int ny, double** ppdView, double** ppdDet,
1090 int iNumDetWithZeros, double dZeropadRatio, double dDetInc)
1092 double dLength = viewDiameter();
1093 // double dLength = phmLen();
1094 double xMin = -dLength / 2;
1095 double xMax = xMin + dLength;
1096 double yMin = -dLength / 2;
1097 double yMax = yMin + dLength;
1098 double xCent = (xMin + xMax) / 2;
1099 double yCent = (yMin + yMax) / 2;
1101 xMin = (xMin - xCent) * dZeropadRatio + xCent;
1102 xMax = (xMax - xCent) * dZeropadRatio + xCent;
1103 yMin = (yMin - yCent) * dZeropadRatio + yCent;
1104 yMax = (yMax - yCent) * dZeropadRatio + yCent;
1106 double xInc = (xMax - xMin) / nx; // size of cells
1107 double yInc = (yMax - yMin) / ny;
1109 // +1 is correct for frequency data, ndet-1 is correct for projections
1110 int iDetCenter = (iNumDetWithZeros - 1) / 2; // index refering to L=0 projection
1111 if (isEven (iNumDetWithZeros))
1112 iDetCenter = (iNumDetWithZeros + 1) / 2;
1114 // Calculates polar coordinates (view#, det#) for each point on phantom grid
1115 double x = xMin + xInc / 2; // Rectang coords of center of pixel
1116 for (unsigned int ix = 0; ix < nx; x += xInc, ix++) {
1117 double y = yMin + yInc / 2;
1118 for (unsigned int iy = 0; iy < ny; y += yInc, iy++) {
1119 double r = ::sqrt (x * x + y * y);
1120 double phi = atan2 (y, x);
1129 ppdView[ix][iy] = (phi - m_rotStart) / m_rotInc;
1130 ppdDet[ix][iy] = (r / dDetInc) + iDetCenter;
1136 Projections::interpolatePolar (ImageFileArray& v, ImageFileArray& vImag,
1137 unsigned int nx, unsigned int ny, std::complex<double>** ppcDetValue, double** ppdView,
1138 double** ppdDet, unsigned int nView, unsigned int nDet, unsigned int nDetWithZeros, int iInterpolationID)
1140 typedef std::complex<double> complexValue;
1142 BilinearInterpolator<complexValue>* pBilinear;
1143 if (iInterpolationID == POLAR_INTERP_BILINEAR)
1144 pBilinear = new BilinearInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
1146 BicubicPolyInterpolator<complexValue>* pBicubic;
1147 if (iInterpolationID == POLAR_INTERP_BICUBIC)
1148 pBicubic = new BicubicPolyInterpolator<complexValue> (ppcDetValue, nView, nDetWithZeros);
1150 for (unsigned int ix = 0; ix < ny; ix++) {
1151 for (unsigned int iy = 0; iy < ny; iy++) {
1153 if (iInterpolationID == POLAR_INTERP_NEAREST) {
1154 unsigned int iView = nearest<int> (ppdView[ix][iy]);
1155 unsigned int iDet = nearest<int> (ppdDet[ix][iy]);
1156 if (iView == nView) {
1158 iDet = m_nDet - iDet;
1160 if (iDet >= 0 && iDet < nDetWithZeros && iView >= 0 && iView < nView) {
1161 v[ix][iy] = ppcDetValue[iView][iDet].real();
1163 vImag[ix][iy] = ppcDetValue[iView][iDet].imag();
1167 } else if (iInterpolationID == POLAR_INTERP_BILINEAR) {
1168 std::complex<double> vInterp = pBilinear->interpolate (ppdView[ix][iy], ppdDet[ix][iy]);
1169 v[ix][iy] = vInterp.real();
1171 vImag[ix][iy] = vInterp.imag();
1172 } else if (iInterpolationID == POLAR_INTERP_BICUBIC) {
1173 std::complex<double> vInterp = pBicubic->interpolate (ppdView[ix][iy], ppdDet[ix][iy]);
1174 v[ix][iy] = vInterp.real();
1176 vImag[ix][iy] = vInterp.imag();
1183 Projections::initFromSomatomAR_STAR (int iNViews, int iNDets, unsigned char* pData, unsigned long lDataLength)
1185 init (iNViews, iNDets);
1186 m_geometry = Scanner::GEOMETRY_EQUIANGULAR;
1187 m_dFocalLength = 510;
1188 m_dSourceDetectorLength = 890;
1189 m_detInc = convertDegreesToRadians (3.06976 / 60);
1190 m_dFanBeamAngle = iNDets * m_detInc;
1191 m_detStart = -(m_dFanBeamAngle / 2);
1192 m_rotInc = TWOPI / static_cast<double>(iNViews);
1194 m_dViewDiameter = sin (m_dFanBeamAngle / 2) * m_dFocalLength * 2;
1196 if (! ((iNViews == 750 && lDataLength == 1560000L) || (iNViews == 950 && lDataLength == 1976000L)
1197 || (iNViews == 1500 && lDataLength == 3120000)))
1200 double dCenter = (iNDets - 1.) / 2.; // change from (Nm+1)/2 because of 0 vs. 1 indexing
1201 double* pdCosScale = new double [iNDets];
1202 for (int i = 0; i < iNDets; i++)
1203 pdCosScale[i] = 1. / cos ((i - dCenter) * m_detInc);
1206 for (int iv = 0; iv < iNViews; iv++) {
1207 unsigned char* pArgBase = pData + lDataPos;
1208 unsigned char* p = pArgBase+0; SwapBytes4IfLittleEndian (p);
1209 long lProjNumber = *reinterpret_cast<long*>(p);
1211 p = pArgBase+20; SwapBytes4IfLittleEndian (p);
1212 long lEscale = *reinterpret_cast<long*>(p);
1214 p = pArgBase+28; SwapBytes4IfLittleEndian (p);
1215 long lTime = *reinterpret_cast<long*>(p);
1217 p = pArgBase + 4; SwapBytes4IfLittleEndian (p);
1218 double dAlpha = *reinterpret_cast<float*>(p) + HALFPI;
1220 p = pArgBase+12; SwapBytes4IfLittleEndian (p);
1221 double dAlign = *reinterpret_cast<float*>(p);
1223 p = pArgBase + 16; SwapBytes4IfLittleEndian (p);
1224 double dMaxValue = *reinterpret_cast<float*>(p);
1226 DetectorArray& detArray = getDetectorArray (iv);
1227 detArray.setViewAngle (dAlpha);
1228 DetectorValue* detval = detArray.detValues();
1230 double dViewScale = 1. / (2294.4871 * ::pow (2.0, -lEscale));
1232 for (int id = 0; id < iNDets; id++) {
1233 int iV = pData[lDataPos+1] + (pData[lDataPos] << 8);
1234 if (iV > 32767) // two's complement signed conversion
1236 detval[id] = iV * dViewScale * pdCosScale[id];
1240 for (int k = iNDets - 2; k >= 0; k--)
1241 detval[k+1] = detval[k];
1251 Projections::interpolateToParallel () const
1253 if (m_geometry == Scanner::GEOMETRY_PARALLEL)
1254 return const_cast<Projections*>(this);
1257 int nView = m_nView;
1258 Projections* pProjNew = new Projections (nView, nDet);
1259 pProjNew->m_geometry = Scanner::GEOMETRY_PARALLEL;
1260 pProjNew->m_dFocalLength = m_dFocalLength;
1261 pProjNew->m_dSourceDetectorLength = m_dSourceDetectorLength;
1262 pProjNew->m_dViewDiameter = m_dViewDiameter;
1263 pProjNew->m_dFanBeamAngle = m_dFanBeamAngle;
1264 pProjNew->m_calcTime = 0;
1265 pProjNew->m_remark = m_remark;
1266 pProjNew->m_remark += "; Interpolate to Parallel";
1267 pProjNew->m_label.setLabelType (Array2dFileLabel::L_HISTORY);
1268 pProjNew->m_label.setLabelString (pProjNew->m_remark);
1269 pProjNew->m_label.setCalcTime (pProjNew->m_calcTime);
1270 pProjNew->m_label.setDateTime (pProjNew->m_year, pProjNew->m_month, pProjNew->m_day, pProjNew->m_hour, pProjNew->m_minute, pProjNew->m_second);
1272 pProjNew->m_rotStart = 0;
1273 #ifdef CONVERT_PARALLEL_PI
1274 pProjNew->m_rotInc = PI / nView;;
1276 pProjNew->m_rotInc = TWOPI / nView;
1278 pProjNew->m_detStart = -m_dViewDiameter / 2;
1279 pProjNew->m_detInc = m_dViewDiameter / nDet;
1280 if (isEven (nDet)) // even
1281 pProjNew->m_detInc = m_dViewDiameter / (nDet - 1);
1283 ParallelRaysums parallel (this, ParallelRaysums::THETA_RANGE_NORMALIZE_TO_TWOPI);
1285 double* pdThetaValuesForT = new double [pProjNew->nView()];
1286 double* pdRaysumsForT = new double [pProjNew->nView()];
1288 // interpolate to evenly spaced theta (views)
1289 double dDetPos = pProjNew->m_detStart;
1290 for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc) {
1291 parallel.getThetaAndRaysumsForT (iD, pdThetaValuesForT, pdRaysumsForT);
1293 double dViewAngle = m_rotStart;
1294 int iLastFloor = -1;
1295 for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
1296 DetectorValue* detValues = pProjNew->getDetectorArray (iV).detValues();
1297 LinearInterpolator<double> interp (pdThetaValuesForT, pdRaysumsForT, pProjNew->nView(), false);
1298 detValues[iD] = interp.interpolate (dViewAngle, &iLastFloor);
1301 delete pdThetaValuesForT;
1302 delete pdRaysumsForT;
1304 // interpolate to evenly space t (detectors)
1305 double* pdOriginalDetPositions = new double [pProjNew->nDet()];
1306 parallel.getDetPositions (pdOriginalDetPositions);
1308 double* pdDetValueCopy = new double [pProjNew->nDet()];
1309 double dViewAngle = m_rotStart;
1310 for (int iV = 0; iV < pProjNew->nView(); iV++, dViewAngle += pProjNew->m_rotInc) {
1311 DetectorArray& detArray = pProjNew->getDetectorArray (iV);
1312 DetectorValue* detValues = detArray.detValues();
1313 detArray.setViewAngle (dViewAngle);
1315 for (int i = 0; i < pProjNew->nDet(); i++)
1316 pdDetValueCopy[i] = detValues[i];
1318 double dDetPos = pProjNew->m_detStart;
1319 int iLastFloor = -1;
1320 LinearInterpolator<double> interp (pdOriginalDetPositions, pdDetValueCopy, pProjNew->nDet(), false);
1321 for (int iD = 0; iD < pProjNew->nDet(); iD++, dDetPos += pProjNew->m_detInc)
1322 detValues[iD] = interp.interpolate (dDetPos, &iLastFloor);
1324 delete pdDetValueCopy;
1325 delete pdOriginalDetPositions;
1331 ///////////////////////////////////////////////////////////////////////////////
1333 // Class ParallelRaysums
1335 // Used for converting divergent beam raysums into Parallel raysums
1337 ///////////////////////////////////////////////////////////////////////////////
1339 ParallelRaysums::ParallelRaysums (const Projections* pProjections, int iThetaRange)
1340 : m_iNumCoordinates(0), m_iNumView(pProjections->nView()), m_iNumDet(pProjections->nDet()),
1341 m_iThetaRange (iThetaRange), m_pCoordinates(NULL)
1343 int iGeometry = pProjections->geometry();
1344 double dDetInc = pProjections->detInc();
1345 double dDetStart = pProjections->detStart();
1346 double dFocalLength = pProjections->focalLength();
1348 m_iNumCoordinates = m_iNumView * m_iNumDet;
1349 m_pCoordinates = new ParallelRaysumCoordinate [m_iNumCoordinates];
1350 m_vecpCoordinates.reserve (m_iNumCoordinates);
1351 for (int i = 0; i < m_iNumCoordinates; i++)
1352 m_vecpCoordinates[i] = m_pCoordinates + i;
1354 int iCoordinate = 0;
1355 for (int iV = 0; iV < m_iNumView; iV++) {
1356 double dViewAngle = pProjections->getDetectorArray(iV).viewAngle();
1357 const DetectorValue* detValues = pProjections->getDetectorArray(iV).detValues();
1359 double dDetPos = dDetStart;
1360 for (int iD = 0; iD < m_iNumDet; iD++) {
1361 ParallelRaysumCoordinate* pC = m_vecpCoordinates[iCoordinate++];
1363 if (iGeometry == Scanner::GEOMETRY_PARALLEL) {
1364 pC->m_dTheta = dViewAngle;
1366 } else if (iGeometry == Scanner::GEOMETRY_EQUILINEAR) {
1367 double dFanAngle = atan (dDetPos / pProjections->sourceDetectorLength());
1368 pC->m_dTheta = dViewAngle + dFanAngle;
1369 pC->m_dT = dFocalLength * sin(dFanAngle);
1371 } else if (iGeometry == Scanner::GEOMETRY_EQUIANGULAR) {
1372 // fan angle is same as dDetPos
1373 pC->m_dTheta = dViewAngle + dDetPos;
1374 pC->m_dT = dFocalLength * sin (dDetPos);
1376 if (m_iThetaRange != THETA_RANGE_UNCONSTRAINED) {
1377 pC->m_dTheta = normalizeAngle (pC->m_dTheta);
1378 if (m_iThetaRange == THETA_RANGE_FOLD_TO_PI && pC->m_dTheta >= PI) {
1380 pC->m_dT = -pC->m_dT;
1383 pC->m_dRaysum = detValues[iD];
1389 ParallelRaysums::~ParallelRaysums()
1391 delete m_pCoordinates;
1394 ParallelRaysums::CoordinateContainer&
1395 ParallelRaysums::getSortedByTheta()
1397 if (m_vecpSortedByTheta.size() == 0) {
1398 m_vecpSortedByTheta.resize (m_iNumCoordinates);
1399 for (int i = 0; i < m_iNumCoordinates; i++)
1400 m_vecpSortedByTheta[i] = m_vecpCoordinates[i];
1401 std::sort (m_vecpSortedByTheta.begin(), m_vecpSortedByTheta.end(), ParallelRaysumCoordinate::compareByTheta);
1404 return m_vecpSortedByTheta;
1407 ParallelRaysums::CoordinateContainer&
1408 ParallelRaysums::getSortedByT()
1410 if (m_vecpSortedByT.size() == 0) {
1411 m_vecpSortedByT.resize (m_iNumCoordinates);
1412 for (int i = 0; i < m_iNumCoordinates; i++)
1413 m_vecpSortedByT[i] = m_vecpCoordinates[i];
1414 std::sort (m_vecpSortedByT.begin(), m_vecpSortedByT.end(), ParallelRaysumCoordinate::compareByT);
1417 return m_vecpSortedByT;
1422 ParallelRaysums::getLimits (double* dMinT, double* dMaxT, double* dMinTheta, double* dMaxTheta) const
1424 if (m_iNumCoordinates <= 0)
1427 *dMinT = *dMaxT = m_vecpCoordinates[0]->m_dT;
1428 *dMinTheta = *dMaxTheta = m_vecpCoordinates[0]->m_dTheta;
1430 for (int i = 0; i < m_iNumCoordinates; i++) {
1431 double dT = m_vecpCoordinates[i]->m_dT;
1432 double dTheta = m_vecpCoordinates[i]->m_dTheta;
1436 else if (dT > *dMaxT)
1439 if (dTheta < *dMinTheta)
1440 *dMinTheta = dTheta;
1441 else if (dTheta > *dMaxTheta)
1442 *dMaxTheta = dTheta;
1447 ParallelRaysums::getThetaAndRaysumsForT (int iTheta, double* pTheta, double* pRaysum)
1449 const CoordinateContainer& coordsT = getSortedByT();
1451 int iBase = iTheta * m_iNumView;
1452 for (int i = 0; i < m_iNumView; i++) {
1453 int iPos = iBase + i;
1454 pTheta[i] = coordsT[iPos]->m_dTheta;
1455 pRaysum[i] = coordsT[iPos]->m_dRaysum;
1460 ParallelRaysums::getDetPositions (double* pdDetPos)
1462 const CoordinateContainer& coordsT = getSortedByT();
1465 for (int i = 0; i < m_iNumDet; i++) {
1466 pdDetPos[i] = coordsT[iPos]->m_dT;