10 #define COLUMNDATA_H 1
34 class ColumnData :
public Column
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 =
"");
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);
57 void data (
int i, T value);
65 ColumnData< T > & operator=(
const ColumnData< T > &right);
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);
73 virtual void insertRows (
long first,
long number = 1);
74 virtual void deleteRows (
long first,
long number = 1);
76 virtual size_t getStoredDataSize()
const;
88 std::vector<T> m_data;
97 inline void ColumnData<T>::readData (
long firstRow,
long nelements,
long firstElem)
99 readColumnData(firstRow,nelements,
static_cast<T*
>(0));
102 template <
typename T>
103 inline const T ColumnData<T>::minLegalValue ()
const
105 return m_minLegalValue;
108 template <
typename T>
109 inline void ColumnData<T>::minLegalValue (T value)
111 m_minLegalValue = value;
114 template <
typename T>
115 inline const T ColumnData<T>::maxLegalValue ()
const
117 return m_maxLegalValue;
120 template <
typename T>
121 inline void ColumnData<T>::maxLegalValue (T value)
123 m_maxLegalValue = value;
126 template <
typename T>
127 inline const T ColumnData<T>::minDataValue ()
const
129 return m_minDataValue;
132 template <
typename T>
133 inline void ColumnData<T>::minDataValue (T value)
135 m_minDataValue = value;
138 template <
typename T>
139 inline const T ColumnData<T>::maxDataValue ()
const
141 return m_maxDataValue;
144 template <
typename T>
145 inline void ColumnData<T>::maxDataValue (T value)
147 m_maxDataValue = value;
150 template <
typename T>
151 inline const std::vector<T>& ColumnData<T>::data ()
const
156 template <
typename T>
157 inline void ColumnData<T>::setData (
const std::vector<T>& value)
162 template <
typename T>
163 inline T ColumnData<T>::data (
int i)
166 return m_data[i - 1];
169 template <
typename T>
170 inline void ColumnData<T>::data (
int i, T value)
173 m_data[i - 1] = value;
178 template <
typename T>
179 ColumnData<T>::ColumnData(
const ColumnData<T> &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),
189 template <
typename T>
190 ColumnData<T>::ColumnData (Table* p)
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),
212 template <
typename T>
213 ColumnData<T>::~ColumnData()
218 template <
typename T>
219 void ColumnData<T>::readColumnData (
long firstRow,
long nelements, T* nullValue)
221 if ( rows() < nelements )
223 std::cerr <<
"CCfits: More data requested than contained in table. ";
224 std::cerr <<
"Extracting complete column.\n";
231 FITSUtil::auto_array_ptr<T> array(
new T[nelements]);
235 if ( fits_read_col(fitsPointer(),type(), index(), firstRow, 1,
236 nelements, nullValue, array.get(), &anynul, &status) )
throw FitsError(status);
239 if (m_data.size() !=
static_cast<size_t>( rows() ) ) m_data.resize(rows());
241 std::copy(&array[0],&array[nelements],m_data.begin()+firstRow-1);
242 if (nelements == rows()) isRead(
true);
245 template <
typename T>
246 bool ColumnData<T>::compare (
const Column &right)
const
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++)
254 if (m_data[i] != that.m_data[i])
return false;
259 template <
typename T>
260 ColumnData<T>* ColumnData<T>::clone ()
const
262 return new ColumnData<T>(*
this);
265 template <
typename T>
266 std::ostream& ColumnData<T>::put (std::ostream& s)
const
269 if (FITS::verboseMode() && type() != Tstring)
271 s <<
" Column Legal limits: ( " << m_minLegalValue <<
"," << m_maxLegalValue <<
" )\n"
272 <<
" Column Data limits: ( " << m_minDataValue <<
"," << m_maxDataValue <<
" )\n";
276 std::ostream_iterator<T> output(s,
"\n");
279 std::copy(m_data.begin(),m_data.end(),output);
285 template <
typename T>
286 void ColumnData<T>::writeData (T* indata,
long nRows,
long firstRow, T* nullValue)
299 if (fits_write_colnull(fitsPointer(), type(), index(), firstRow, 1, nRows,
300 indata, nullValue, &status) != 0)
throw FitsError(status);
304 if (fits_write_col(fitsPointer(), type(), index(), firstRow, 1, nRows,
305 indata, &status) != 0)
throw FitsError(status);
311 if (status == NO_NULL)
throw NoNullValue(name());
314 long elementsToWrite(nRows + firstRow -1);
316 if (elementsToWrite >
static_cast<long>(m_data.size()))
318 m_data.resize(elementsToWrite,T());
321 std::copy(&indata[0],&indata[nRows],m_data.begin()+firstRow-1);
323 parent()->updateRows();
327 template <
typename T>
328 void ColumnData<T>::writeData (
const std::vector<T>& indata,
long firstRow, T* nullValue)
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);
336 template <
typename T>
337 void ColumnData<T>::insertRows (
long first,
long number)
345 if (first >= 0 && first <=
static_cast<long>(m_data.size()))
347 typename std::vector<T>::iterator in;
350 in = m_data.begin()+first;
358 m_data.insert(in,number,T());
362 template <
typename T>
363 void ColumnData<T>::deleteRows (
long first,
long number)
369 const long curSize =
static_cast<long>(m_data.size());
370 if (curSize>0 && first <= curSize)
372 const long last = std::min(curSize, first-1+number);
373 m_data.erase(m_data.begin()+first-1,m_data.begin()+last);
377 template <
typename T>
378 size_t ColumnData<T>::getStoredDataSize()
const
380 return m_data.size();
383 template <
typename T>
384 void ColumnData<T>::setDataLimits (T* limits)
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]);
397 #if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
399 inline void ColumnData<complex<float> >::setDataLimits (complex<float>* limits)
401 m_minLegalValue = limits[0];
402 m_maxLegalValue = limits[1];
403 m_minDataValue = limits[2];
404 m_maxDataValue = limits[3];
408 void ColumnData<complex<float> >::setDataLimits (complex<float>* limits);
411 #if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
413 inline void ColumnData<complex<double> >::setDataLimits (complex<double>* limits)
415 m_minLegalValue = limits[0];
416 m_maxLegalValue = limits[1];
417 m_minDataValue = limits[2];
418 m_maxDataValue = limits[3];
422 void ColumnData<complex<double> >::setDataLimits (complex<double>* limits);
426 #if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
428 inline void ColumnData<string>::readColumnData (
long firstRow,
434 throw Column::InvalidNumberOfRows((
int)nelements);
435 if (firstRow < 1 || (firstRow+nelements-1)>rows())
436 throw Column::InvalidRowNumber(name());
441 char** array =
new char*[nelements];
444 for (
long i=0; i<nelements; ++i)
445 array[i]=
static_cast<char*
>(0);
446 bool isError =
false;
453 nulval =
const_cast<char*
>(nullValue->c_str());
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))
474 for (
long j=0; j<nelements; ++j)
476 array[j] =
new char[strLengths[j] + 1];
479 const long lastRow = firstRow+nelements-1;
480 for (
long iRow=firstRow; !isError && iRow<=lastRow; ++iRow)
482 if (fits_read_col_str(fitsPointer(),index(), iRow, 1, 1,
483 nulval, &array[iRow-firstRow], &anynul,&status) )
487 delete [] strLengths;
493 for (
long j=0; j<nelements; ++j)
495 array[j] =
new char[width() + 1];
497 if (fits_read_col_str(fitsPointer(),index(), firstRow,1,nelements,
498 nulval,array, &anynul,&status))
507 for (
long j = 0; j < nelements; ++j)
513 throw FitsError(status);
516 if (m_data.size() !=
static_cast<size_t>(rows()))
517 setData(std::vector<String>(rows(),String(nulval)));
519 for (
long j=0; j<nelements; ++j)
521 m_data[j - 1 + firstRow] = String(array[j]);
524 for (
long j=0; j<nelements; j++)
530 if (nelements == rows()) isRead(
true);
535 void ColumnData<string>::readColumnData (
long firstRow,
long nelements,
string* nullValue);
539 #if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
541 inline void ColumnData<complex<float> >::readColumnData (
long firstRow,
543 complex<float>* nullValue)
548 FITSUtil::auto_array_ptr<float> pArray(
new float[nelements*2]);
549 float* array = pArray.get();
554 if (fits_read_col_cmp(fitsPointer(),index(), firstRow,1,nelements,
555 nulval,array, &anynul,&status) )
throw FitsError(status);
558 if (m_data.size() != rows()) m_data.resize(rows());
562 for (
int j = 0; j < nelements; ++j)
565 m_data[j - 1 + firstRow] = std::complex<float>(array[2*j],array[2*j+1]);
567 if (nelements == rows()) isRead(
true);
572 void ColumnData<complex<float> >::readColumnData (
long firstRow,
long nelements,complex<float>* nullValue );
575 #if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
577 inline void ColumnData<complex<double> >::readColumnData (
long firstRow,
579 complex<double>* nullValue)
584 FITSUtil::auto_array_ptr<double> pArray(
new double[nelements*2]);
585 double* array = pArray.get();
590 if (fits_read_col_dblcmp(fitsPointer(), index(), firstRow,1,nelements,
591 nulval,array, &anynul,&status) )
throw FitsError(status);
596 if (m_data.size() != rows()) setData(std::vector<complex<double> >(rows(),nulval));
600 for (
int j = 0; j < nelements; j++)
603 m_data[j - 1 + firstRow] = std::complex<double>(array[2*j],array[2*j+1]);
605 if (nelements == rows()) isRead(
true);
610 void ColumnData<complex<double> >::readColumnData (
long firstRow,
long nelements,complex<double>* nullValue);
613 #if SPEC_TEMPLATE_DECL_DEFECT
615 inline void ColumnData<string>::writeData (
const std::vector<string>& indata,
616 long firstRow,
string* nullValue)
619 char** columnData=FITSUtil::CharArray(indata);
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)
628 m_data.resize(elementsToWrite,
"");
629 std::copy(__tmp.begin(),__tmp.end(),m_data.begin());
631 std::copy(indata.begin(),indata.end(),m_data.begin()+firstRow-1);
634 for (
size_t i = 0; i < indata.size(); ++i)
636 delete [] columnData[i];
638 delete [] columnData;
642 void ColumnData<string>::writeData (
const std::vector<string>& inData,
long firstRow,
string* nullValue);
645 #ifdef SPEC_TEMPLATE_DECL_DEFECT
647 inline void ColumnData<complex<float> >::writeData (
const std::vector<complex<float> >& inData,
649 complex<float>* nullValue)
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)
658 Data[ 2*j] = inData[j].real();
659 Data[ 2*j + 1] = inData[j].imag();
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()))
671 m_data.resize(elementsToWrite);
674 std::copy(inData.begin(),inData.end(),m_data.begin()+firstRow-1);
677 parent()->updateRows();
682 m_data.resize(__tmp.size());
690 void ColumnData<complex<float> >::writeData (
const std::vector<complex<float> >& inData,
long firstRow,
691 complex<float>* nullValue);
694 #ifdef SPEC_TEMPLATE_DECL_DEFECT
696 inline void ColumnData<complex<double> >::writeData (
const std::vector<complex<double> >& inData,
698 complex<double>* nullValue)
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)
707 pData[ 2*j] = inData[j].real();
708 pData[ 2*j + 1] = inData[j].imag();
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()))
720 m_data.resize(elementsToWrite);
723 std::copy(inData.begin(),inData.end(),m_data.begin()+firstRow-1);
726 parent()->updateRows();
731 m_data.resize(__tmp.size());
739 void ColumnData<complex<double> >::writeData (
const std::vector<complex<double> >& inData,
long firstRow,
740 complex<double>* nullValue);
Column(const Column &right)
copy constructor, used in copying Columns to standard library containers.
Definition: Column.cxx:171
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