CCfits 2.7
ColumnData.h
1// Astrophysics Science Division,
2// NASA/ Goddard Space Flight Center
3// HEASARC
4// http://heasarc.gsfc.nasa.gov
5// e-mail: ccfits@legacy.gsfc.nasa.gov
6//
7// Original author: Ben Dorman
8
9#ifndef COLUMNDATA_H
10#define COLUMNDATA_H 1
11#include "CCfits.h"
12
13// vector
14#include <vector>
15// Column
16#include "Column.h"
17#ifdef _MSC_VER
18#include "MSconfig.h"
19#endif
20
21#include <complex>
22#include <memory>
23#include <iterator>
24#include "FITSUtil.h"
25using std::complex;
26#include "FITS.h"
27
28
29namespace CCfits {
30
31
32
33 template <typename T>
34 class ColumnData : public Column //## Inherits: <unnamed>%385E51565EE8
35 {
36
37 public:
38 ColumnData(const ColumnData< T > &right);
39 ColumnData (Table* p = 0);
40 ColumnData (int columnIndex, const string &columnName, ValueType type, const String &format, const String &unit, Table* p, int rpt = 1, long w = 1, const String &comment = "");
41 ~ColumnData();
42
43 virtual ColumnData<T>* clone () const;
44 virtual void readData (long firstRow, long nelements, long firstElem = 1);
45 void setDataLimits (T* limits);
46 const T minLegalValue () const;
47 void minLegalValue (T value);
48 const T maxLegalValue () const;
49 void maxLegalValue (T value);
50 const T minDataValue () const;
51 void minDataValue (T value);
52 const T maxDataValue () const;
53 void maxDataValue (T value);
54 const std::vector<T>& data () const;
55 void setData (const std::vector<T>& value);
56 T data (int i);
57 void data (int i, T value);
58
59 // Additional Public Declarations
60 friend class Column;
61 protected:
62 // Additional Protected Declarations
63
64 private:
65 ColumnData< T > & operator=(const ColumnData< T > &right);
66
67 void readColumnData (long firstRow, long nelements, T* nullValue = 0);
68 virtual bool compare (const Column &right) const;
69 virtual std::ostream& put (std::ostream& s) const;
70 void writeData (T* indata, long nRows = 1, long firstRow = 1, T* nullValue = 0);
71 void writeData (const std::vector<T>& indata, long firstRow = 1, T* nullValue = 0);
72 // Insert one or more blank rows into a FITS column.
73 virtual void insertRows (long first, long number = 1);
74 virtual void deleteRows (long first, long number = 1);
75
76 virtual size_t getStoredDataSize() const;
77
78 // Additional Private Declarations
79
80 private: //## implementation
81 // Data Members for Class Attributes
82 T m_minLegalValue;
83 T m_maxLegalValue;
84 T m_minDataValue;
85 T m_maxDataValue;
86
87 // Data Members for Associations
88 std::vector<T> m_data;
89
90 // Additional Implementation Declarations
91
92 };
93
94 // Parameterized Class CCfits::ColumnData
95
96 template <typename T>
97 inline void ColumnData<T>::readData (long firstRow, long nelements, long firstElem)
98 {
99 readColumnData(firstRow,nelements,static_cast<T*>(0));
100 }
101
102 template <typename T>
103 inline const T ColumnData<T>::minLegalValue () const
104 {
105 return m_minLegalValue;
106 }
107
108 template <typename T>
109 inline void ColumnData<T>::minLegalValue (T value)
110 {
111 m_minLegalValue = value;
112 }
113
114 template <typename T>
115 inline const T ColumnData<T>::maxLegalValue () const
116 {
117 return m_maxLegalValue;
118 }
119
120 template <typename T>
121 inline void ColumnData<T>::maxLegalValue (T value)
122 {
123 m_maxLegalValue = value;
124 }
125
126 template <typename T>
127 inline const T ColumnData<T>::minDataValue () const
128 {
129 return m_minDataValue;
130 }
131
132 template <typename T>
133 inline void ColumnData<T>::minDataValue (T value)
134 {
135 m_minDataValue = value;
136 }
137
138 template <typename T>
139 inline const T ColumnData<T>::maxDataValue () const
140 {
141 return m_maxDataValue;
142 }
143
144 template <typename T>
145 inline void ColumnData<T>::maxDataValue (T value)
146 {
147 m_maxDataValue = value;
148 }
149
150 template <typename T>
151 inline const std::vector<T>& ColumnData<T>::data () const
152 {
153 return m_data;
154 }
155
156 template <typename T>
157 inline void ColumnData<T>::setData (const std::vector<T>& value)
158 {
159 m_data = value;
160 }
161
162 template <typename T>
163 inline T ColumnData<T>::data (int i)
164 {
165 // return data stored in the ith row, which is in the i-1 th location in the array.
166 return m_data[i - 1];
167 }
168
169 template <typename T>
170 inline void ColumnData<T>::data (int i, T value)
171 {
172 // assign data to i-1 th location in the array, representing the ith row.
173 m_data[i - 1] = value;
174 }
175
176 // Parameterized Class CCfits::ColumnData
177
178 template <typename T>
179 ColumnData<T>::ColumnData(const ColumnData<T> &right)
180 :Column(right),
181 m_minLegalValue(right.m_minLegalValue),
182 m_maxLegalValue(right.m_maxLegalValue),
183 m_minDataValue(right.m_minDataValue),
184 m_maxDataValue(right.m_maxDataValue),
185 m_data(right.m_data)
186 {
187 }
188
189 template <typename T>
190 ColumnData<T>::ColumnData (Table* p)
191 : Column(p),
192 m_minLegalValue(),
193 m_maxLegalValue(),
194 m_minDataValue(),
195 m_maxDataValue(),
196 m_data()
197 {
198 }
199
200 template <typename T>
201 ColumnData<T>::ColumnData (int columnIndex, const string &columnName, ValueType type, const String &format, const String &unit, Table* p, int rpt, long w, const String &comment)
202 : Column(columnIndex,columnName,type,format,unit,p,rpt,w,comment),
203 m_minLegalValue(),
204 m_maxLegalValue(),
205 m_minDataValue(),
206 m_maxDataValue(),
207 m_data()
208 {
209 }
210
211
212 template <typename T>
213 ColumnData<T>::~ColumnData()
214 {
215 }
216
217
218 template <typename T>
219 void ColumnData<T>::readColumnData (long firstRow, long nelements, T* nullValue)
220 {
221 if ( rows() < nelements )
222 {
223 std::cerr << "CCfits: More data requested than contained in table. ";
224 std::cerr << "Extracting complete column.\n";
225 nelements = rows();
226 }
227
228 int status(0);
229 int anynul(0);
230
231 FITSUtil::auto_array_ptr<T> array(new T[nelements]);
232
233 makeHDUCurrent();
234
235 if ( fits_read_col(fitsPointer(),type(), index(), firstRow, 1,
236 nelements, nullValue, array.get(), &anynul, &status) ) throw FitsError(status);
237
238
239 if (m_data.size() != static_cast<size_t>( rows() ) ) m_data.resize(rows());
240
241 std::copy(&array[0],&array[nelements],m_data.begin()+firstRow-1);
242 if (nelements == rows()) isRead(true);
243 }
244
245 template <typename T>
246 bool ColumnData<T>::compare (const Column &right) const
247 {
248 if ( !Column::compare(right) ) return false;
249 const ColumnData<T>& that = static_cast<const ColumnData<T>&>(right);
250 unsigned int n = m_data.size();
251 if ( that.m_data.size() != n ) return false;
252 for (unsigned int i = 0; i < n ; i++)
253 {
254 if (m_data[i] != that.m_data[i]) return false;
255 }
256 return true;
257 }
258
259 template <typename T>
260 ColumnData<T>* ColumnData<T>::clone () const
261 {
262 return new ColumnData<T>(*this);
263 }
264
265 template <typename T>
266 std::ostream& ColumnData<T>::put (std::ostream& s) const
267 {
268 Column::put(s);
269 if (FITS::verboseMode() && type() != Tstring)
270 {
271 s << " Column Legal limits: ( " << m_minLegalValue << "," << m_maxLegalValue << " )\n"
272 << " Column Data limits: ( " << m_minDataValue << "," << m_maxDataValue << " )\n";
273 }
274 if (!m_data.empty())
275 {
276 std::ostream_iterator<T> output(s,"\n");
277 // output each row on a separate line.
278 // user can supply manipulators to stream for formatting.
279 std::copy(m_data.begin(),m_data.end(),output);
280 }
281
282 return s;
283 }
284
285 template <typename T>
286 void ColumnData<T>::writeData (T* indata, long nRows, long firstRow, T* nullValue)
287 {
288
289 // set columnData's data member to equal what's written to file.
290 // indata has size nRows: elements firstRow to firstRow + nRows - 1 will be written.
291 // if this exceeds the current rowlength of the HDU, update the return value for
292 // rows() in the parent after the fitsio call.
293 int status(0);
294
295 try
296 {
297 if (nullValue)
298 {
299 if (fits_write_colnull(fitsPointer(), type(), index(), firstRow, 1, nRows,
300 indata, nullValue, &status) != 0) throw FitsError(status);
301 }
302 else
303 {
304 if (fits_write_col(fitsPointer(), type(), index(), firstRow, 1, nRows,
305 indata, &status) != 0) throw FitsError(status);
306 }
307
308 }
309 catch (FitsError) // the only thing that can throw here.
310 {
311 if (status == NO_NULL) throw NoNullValue(name());
312 else throw;
313 }
314 long elementsToWrite(nRows + firstRow -1);
315
316 if (elementsToWrite > static_cast<long>(m_data.size()))
317 {
318 m_data.resize(elementsToWrite,T());
319 }
320
321 std::copy(&indata[0],&indata[nRows],m_data.begin()+firstRow-1);
322 // tell the Table that the number of rows has changed
323 parent()->updateRows();
324
325 }
326
327 template <typename T>
328 void ColumnData<T>::writeData (const std::vector<T>& indata, long firstRow, T* nullValue)
329 {
330 FITSUtil::CVarray<T> convert;
331 FITSUtil::auto_array_ptr<T> pcolData (convert(indata));
332 T* columnData = pcolData.get();
333 writeData(columnData,indata.size(),firstRow,nullValue);
334 }
335
336 template <typename T>
337 void ColumnData<T>::insertRows (long first, long number)
338 {
339 // Don't assume the calling routine (ie. Table's insertRows)
340 // knows Column's current m_data size. m_data may still be
341 // size 0 if no read operations have been performed on Column.
342 // Therefore perform range checking before inserting.
343
344 // As with CFITSIO, rows are inserted AFTER 1-based 'first'.
345 if (first >= 0 && first <= static_cast<long>(m_data.size()))
346 {
347 typename std::vector<T>::iterator in;
348 if (first !=0)
349 {
350 in = m_data.begin()+first;
351 }
352 else
353 {
354 in = m_data.begin();
355 }
356
357 // non-throwing operations.
358 m_data.insert(in,number,T());
359 }
360 }
361
362 template <typename T>
363 void ColumnData<T>::deleteRows (long first, long number)
364 {
365 // Don't assume the calling routine (ie. Table's deleteRows)
366 // knows Column's current m_data size. m_data may still be
367 // size 0 if no read operations have been performed on Column.
368 // Therefore perform range checking before erasing.
369 const long curSize = static_cast<long>(m_data.size());
370 if (curSize>0 && first <= curSize)
371 {
372 const long last = std::min(curSize, first-1+number);
373 m_data.erase(m_data.begin()+first-1,m_data.begin()+last);
374 }
375 }
376
377 template <typename T>
378 size_t ColumnData<T>::getStoredDataSize() const
379 {
380 return m_data.size();
381 }
382
383 template <typename T>
384 void ColumnData<T>::setDataLimits (T* limits)
385 {
386 m_minLegalValue = limits[0];
387 m_maxLegalValue = limits[1];
388 m_minDataValue = std::max(limits[2],limits[0]);
389 m_maxDataValue = std::min(limits[3],limits[1]);
390 }
391
392 // Additional Declarations
393
394 // all functions that operate on strings or complex data that call cfitsio
395 // need to be specialized.
396
397#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
398template <>
399inline void ColumnData<complex<float> >::setDataLimits (complex<float>* limits)
400 {
401 m_minLegalValue = limits[0];
402 m_maxLegalValue = limits[1];
403 m_minDataValue = limits[2];
404 m_maxDataValue = limits[3];
405 }
406#else
407template <>
408 void ColumnData<complex<float> >::setDataLimits (complex<float>* limits);
409#endif
410
411#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
412template <>
413inline void ColumnData<complex<double> >::setDataLimits (complex<double>* limits)
414 {
415 m_minLegalValue = limits[0];
416 m_maxLegalValue = limits[1];
417 m_minDataValue = limits[2];
418 m_maxDataValue = limits[3];
419 }
420#else
421 template <>
422 void ColumnData<complex<double> >::setDataLimits (complex<double>* limits);
423#endif
424
425
426#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
427 template <>
428 inline void ColumnData<string>::readColumnData (long firstRow,
429 long nelements,
430 string* nullValue)
431 {
432 // nelements = nrows to read.
433 if (nelements < 1)
434 throw Column::InvalidNumberOfRows((int)nelements);
435 if (firstRow < 1 || (firstRow+nelements-1)>rows())
436 throw Column::InvalidRowNumber(name());
437
438 int status = 0;
439 int anynul = 0;
440
441 char** array = new char*[nelements];
442 // Initialize pointers to NULL so we can safely delete
443 // during error handling, even if they haven't been allocated.
444 for (long i=0; i<nelements; ++i)
445 array[i]=static_cast<char*>(0);
446 bool isError = false;
447
448 // Strings are unusual. The variable length case is still
449 // handled by a ColumnData class, not a ColumnVectorData.
450 char* nulval = 0;
451 if (nullValue)
452 {
453 nulval = const_cast<char*>(nullValue->c_str());
454 }
455 else
456 {
457 nulval = new char;
458 *nulval = '\0';
459 }
460 makeHDUCurrent();
461 if (varLength())
462 {
463 long* strLengths = new long[nelements];
464 long* offsets = new long[nelements];
465 if (fits_read_descripts(fitsPointer(), index(), firstRow,
466 nelements, strLengths, offsets, &status))
467 {
468 isError = true;
469 }
470 else
471 {
472 // For variable length cols, must read 1 and only 1 row
473 // at a time into array.
474 for (long j=0; j<nelements; ++j)
475 {
476 array[j] = new char[strLengths[j] + 1];
477 }
478
479 const long lastRow = firstRow+nelements-1;
480 for (long iRow=firstRow; !isError && iRow<=lastRow; ++iRow)
481 {
482 if (fits_read_col_str(fitsPointer(),index(), iRow, 1, 1,
483 nulval, &array[iRow-firstRow], &anynul,&status) )
484 isError=true;
485 }
486 }
487 delete [] strLengths;
488 delete [] offsets;
489 }
490 else
491 {
492 // Fixed length strings, length is stored in Column's m_width.
493 for (long j=0; j<nelements; ++j)
494 {
495 array[j] = new char[width() + 1];
496 }
497 if (fits_read_col_str(fitsPointer(),index(), firstRow,1,nelements,
498 nulval,array, &anynul,&status))
499 isError=true;
500 }
501
502 if (isError)
503 {
504 // It's OK to do this even if error occurred before
505 // array rows were allocated. In that case their pointers
506 // were set to NULL.
507 for (long j = 0; j < nelements; ++j)
508 {
509 delete [] array[j];
510 }
511 delete [] array;
512 delete nulval;
513 throw FitsError(status);
514 }
515
516 if (m_data.size() != static_cast<size_t>(rows()))
517 setData(std::vector<String>(rows(),String(nulval)));
518
519 for (long j=0; j<nelements; ++j)
520 {
521 m_data[j - 1 + firstRow] = String(array[j]);
522 }
523
524 for (long j=0; j<nelements; j++)
525 {
526 delete [] array[j];
527 }
528 delete [] array;
529 delete nulval;
530 if (nelements == rows()) isRead(true);
531
532 }
533#else
534 template <>
535void ColumnData<string>::readColumnData (long firstRow, long nelements, string* nullValue);
536#endif
537
538
539#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
540 template <>
541 inline void ColumnData<complex<float> >::readColumnData (long firstRow,
542 long nelements,
543 complex<float>* nullValue)
544 {
545 // specialization for ColumnData<string>
546 int status(0);
547 int anynul(0);
548 FITSUtil::auto_array_ptr<float> pArray(new float[nelements*2]);
549 float* array = pArray.get();
550 float nulval(0);
551 makeHDUCurrent();
552
553
554 if (fits_read_col_cmp(fitsPointer(),index(), firstRow,1,nelements,
555 nulval,array, &anynul,&status) ) throw FitsError(status);
556
557
558 if (m_data.size() != rows()) m_data.resize(rows());
559
560 // the 'j -1 ' converts to zero based indexing.
561
562 for (int j = 0; j < nelements; ++j)
563 {
564
565 m_data[j - 1 + firstRow] = std::complex<float>(array[2*j],array[2*j+1]);
566 }
567 if (nelements == rows()) isRead(true);
568
569 }
570#else
571template <>
572void ColumnData<complex<float> >::readColumnData (long firstRow, long nelements,complex<float>* nullValue );
573#endif
574
575#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
576 template <>
577 inline void ColumnData<complex<double> >::readColumnData (long firstRow,
578 long nelements,
579 complex<double>* nullValue)
580 {
581 // specialization for ColumnData<complex<double> >
582 int status(0);
583 int anynul(0);
584 FITSUtil::auto_array_ptr<double> pArray(new double[nelements*2]);
585 double* array = pArray.get();
586 double nulval(0);
587 makeHDUCurrent();
588
589
590 if (fits_read_col_dblcmp(fitsPointer(), index(), firstRow,1,nelements,
591 nulval,array, &anynul,&status) ) throw FitsError(status);
592
593
594
595
596 if (m_data.size() != rows()) setData(std::vector<complex<double> >(rows(),nulval));
597
598 // the 'j -1 ' converts to zero based indexing.
599
600 for (int j = 0; j < nelements; j++)
601 {
602
603 m_data[j - 1 + firstRow] = std::complex<double>(array[2*j],array[2*j+1]);
604 }
605 if (nelements == rows()) isRead(true);
606
607 }
608#else
609template <>
610void ColumnData<complex<double> >::readColumnData (long firstRow, long nelements,complex<double>* nullValue);
611#endif
612
613#if SPEC_TEMPLATE_DECL_DEFECT
614 template <>
615 inline void ColumnData<string>::writeData (const std::vector<string>& indata,
616 long firstRow, string* nullValue)
617 {
618 int status=0;
619 char** columnData=FITSUtil::CharArray(indata);
620
621 if ( fits_write_colnull(fitsPointer(), TSTRING, index(), firstRow, 1, indata.size(),
622 columnData, 0, &status) != 0 )
623 throw FitsError(status);
624 unsigned long elementsToWrite (indata.size() + firstRow - 1);
625 std::vector<string> __tmp(m_data);
626 if (m_data.size() < elementsToWrite)
627 {
628 m_data.resize(elementsToWrite,"");
629 std::copy(__tmp.begin(),__tmp.end(),m_data.begin());
630 }
631 std::copy(indata.begin(),indata.end(),m_data.begin()+firstRow-1);
632
633
634 for (size_t i = 0; i < indata.size(); ++i)
635 {
636 delete [] columnData[i];
637 }
638 delete [] columnData;
639 }
640#else
641template <>
642void ColumnData<string>::writeData (const std::vector<string>& inData, long firstRow, string* nullValue);
643#endif
644
645#ifdef SPEC_TEMPLATE_DECL_DEFECT
646 template <>
647 inline void ColumnData<complex<float> >::writeData (const std::vector<complex<float> >& inData,
648 long firstRow,
649 complex<float>* nullValue)
650 {
651 int status(0);
652 int nRows (inData.size());
653 FITSUtil::auto_array_ptr<float> pData(new float[nRows*2]);
654 float* Data = pData.get();
655 std::vector<complex<float> > __tmp(m_data);
656 for (int j = 0; j < nRows; ++j)
657 {
658 Data[ 2*j] = inData[j].real();
659 Data[ 2*j + 1] = inData[j].imag();
660 }
661
662 try
663 {
664
665 if (fits_write_col_cmp(fitsPointer(), index(), firstRow, 1,
666 nRows,Data, &status) != 0) throw FitsError(status);
667 long elementsToWrite(nRows + firstRow -1);
668 if (elementsToWrite > static_cast<long>(m_data.size()))
669 {
670
671 m_data.resize(elementsToWrite);
672 }
673
674 std::copy(inData.begin(),inData.end(),m_data.begin()+firstRow-1);
675
676 // tell the Table that the number of rows has changed
677 parent()->updateRows();
678 }
679 catch (FitsError) // the only thing that can throw here.
680 {
681 // reset to original content and rethrow the exception.
682 m_data.resize(__tmp.size());
683 m_data = __tmp;
684 }
685
686 }
687
688#else
689template <>
690void ColumnData<complex<float> >::writeData (const std::vector<complex<float> >& inData, long firstRow,
691 complex<float>* nullValue);
692#endif
693
694#ifdef SPEC_TEMPLATE_DECL_DEFECT
695 template <>
696 inline void ColumnData<complex<double> >::writeData (const std::vector<complex<double> >& inData,
697 long firstRow,
698 complex<double>* nullValue)
699 {
700 int status(0);
701 int nRows (inData.size());
702 FITSUtil::auto_array_ptr<double> pData(new double[nRows*2]);
703 double* Data = pData.get();
704 std::vector<complex<double> > __tmp(m_data);
705 for (int j = 0; j < nRows; ++j)
706 {
707 pData[ 2*j] = inData[j].real();
708 pData[ 2*j + 1] = inData[j].imag();
709 }
710
711 try
712 {
713
714 if (fits_write_col_dblcmp(fitsPointer(), index(), firstRow, 1,
715 nRows,Data, &status) != 0) throw FitsError(status);
716 long elementsToWrite(nRows + firstRow -1);
717 if (elementsToWrite > static_cast<long>(m_data.size()))
718 {
719
720 m_data.resize(elementsToWrite);
721 }
722
723 std::copy(inData.begin(),inData.end(),m_data.begin()+firstRow-1);
724
725 // tell the Table that the number of rows has changed
726 parent()->updateRows();
727 }
728 catch (FitsError) // the only thing that can throw here.
729 {
730 // reset to original content and rethrow the exception.
731 m_data.resize(__tmp.size());
732 m_data = __tmp;
733 }
734
735 }
736
737#else
738template <>
739void ColumnData<complex<double> >::writeData (const std::vector<complex<double> >& inData, long firstRow,
740 complex<double>* nullValue);
741
742#endif
743} // namespace CCfits
744
745
746#endif
Column(const Column &right)
copy constructor, used in copying Columns to standard library containers.
Definition Column.cxx:171
const String & unit() const
get units of data in Column (TUNITn keyword)
Definition Column.h:1527
const String & comment() const
retrieve comment for Column
Definition Column.h:1517
const String & format() const
return TFORMn keyword
Definition Column.h:1522
ValueType type() const
returns the data type of the column
Definition Column.h:1437
Namespace enclosing all CCfits classes and globals definitions.
Definition AsciiTable.cxx:26
ValueType
CCfits value types and their CFITSIO equivalents (in caps)
Definition CCfits.h:81