This file is indexed.

/usr/include/CCfits/ColumnData.h is in libccfits-dev 2.5+dfsg-1+b1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
//	Astrophysics Science Division,
//	NASA/ Goddard Space Flight Center
//	HEASARC
//	http://heasarc.gsfc.nasa.gov
//	e-mail: ccfits@legacy.gsfc.nasa.gov
//
//	Original author: Ben Dorman

#ifndef COLUMNDATA_H
#define COLUMNDATA_H 1
#include "CCfits.h"

// vector
#include <vector>
// Column
#include "Column.h"
#ifdef _MSC_VER
#include "MSconfig.h"
#endif

#include <complex>
#include <memory>
#include <iterator>
#include "FITSUtil.h"
using std::complex;
#include "FITS.h"


namespace CCfits {



  template <typename T>
  class ColumnData : public Column  //## Inherits: <unnamed>%385E51565EE8
  {

    public:
        ColumnData(const ColumnData< T > &right);
        ColumnData (Table* p = 0);
        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 = "");
        ~ColumnData();

        virtual ColumnData<T>* clone () const;
        virtual void readData (long firstRow, long nelements, long firstElem = 1);
        void setDataLimits (T* limits);
        const T minLegalValue () const;
        void minLegalValue (T value);
        const T maxLegalValue () const;
        void maxLegalValue (T value);
        const T minDataValue () const;
        void minDataValue (T value);
        const T maxDataValue () const;
        void maxDataValue (T value);
        const std::vector<T>& data () const;
        void setData (const std::vector<T>& value);
        T data (int i);
        void data (int i, T value);

      // Additional Public Declarations
        friend class Column;
    protected:
      // Additional Protected Declarations

    private:
        ColumnData< T > & operator=(const ColumnData< T > &right);

        void readColumnData (long firstRow, long nelements, T* nullValue = 0);
        virtual bool compare (const Column &right) const;
        virtual std::ostream& put (std::ostream& s) const;
        void writeData (T* indata, long nRows = 1, long firstRow = 1, T* nullValue = 0);
        void writeData (const std::vector<T>& indata, long firstRow = 1, T* nullValue = 0);
        //	Insert one or more blank rows into a FITS column.
        virtual void insertRows (long first, long number = 1);
        virtual void deleteRows (long first, long number = 1);

      // Additional Private Declarations

    private: //## implementation
      // Data Members for Class Attributes
        T m_minLegalValue;
        T m_maxLegalValue;
        T m_minDataValue;
        T m_maxDataValue;

      // Data Members for Associations
        std::vector<T> m_data;

      // Additional Implementation Declarations

  };

  // Parameterized Class CCfits::ColumnData 

  template <typename T>
  inline void ColumnData<T>::readData (long firstRow, long nelements, long firstElem)
  {
   readColumnData(firstRow,nelements,static_cast<T*>(0));
  }

  template <typename T>
  inline const T ColumnData<T>::minLegalValue () const
  {
    return m_minLegalValue;
  }

  template <typename T>
  inline void ColumnData<T>::minLegalValue (T value)
  {
    m_minLegalValue = value;
  }

  template <typename T>
  inline const T ColumnData<T>::maxLegalValue () const
  {
    return m_maxLegalValue;
  }

  template <typename T>
  inline void ColumnData<T>::maxLegalValue (T value)
  {
    m_maxLegalValue = value;
  }

  template <typename T>
  inline const T ColumnData<T>::minDataValue () const
  {
    return m_minDataValue;
  }

  template <typename T>
  inline void ColumnData<T>::minDataValue (T value)
  {
    m_minDataValue = value;
  }

  template <typename T>
  inline const T ColumnData<T>::maxDataValue () const
  {
    return m_maxDataValue;
  }

  template <typename T>
  inline void ColumnData<T>::maxDataValue (T value)
  {
    m_maxDataValue = value;
  }

  template <typename T>
  inline const std::vector<T>& ColumnData<T>::data () const
  {
    return m_data;
  }

  template <typename T>
  inline void ColumnData<T>::setData (const std::vector<T>& value)
  {
    m_data = value;
  }

  template <typename T>
  inline T ColumnData<T>::data (int i)
  {
    // return data stored in the ith row, which is in the i-1 th location in the array.
    return m_data[i - 1];
  }

  template <typename T>
  inline void ColumnData<T>::data (int i, T value)
  {
    // assign data to i-1 th location in the array, representing the ith row.
    m_data[i - 1] = value;
  }

  // Parameterized Class CCfits::ColumnData 

  template <typename T>
  ColumnData<T>::ColumnData(const ColumnData<T> &right)
      :Column(right),
       m_minLegalValue(right.m_minLegalValue),
       m_maxLegalValue(right.m_maxLegalValue),
       m_minDataValue(right.m_minDataValue),
       m_maxDataValue(right.m_maxDataValue),
       m_data(right.m_data)
  {
  }

  template <typename T>
  ColumnData<T>::ColumnData (Table* p)
  : Column(p),
       m_minLegalValue(),
       m_maxLegalValue(),
       m_minDataValue(),
       m_maxDataValue(), 
       m_data()
  {
  }

  template <typename T>
  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)
        : Column(columnIndex,columnName,type,format,unit,p,rpt,w,comment), 
        m_minLegalValue(),
        m_maxLegalValue(),
        m_minDataValue(),
        m_maxDataValue(),
        m_data()
  {
  }


  template <typename T>
  ColumnData<T>::~ColumnData()
  {
  }


  template <typename T>
  void ColumnData<T>::readColumnData (long firstRow, long nelements, T* nullValue)
  {
  if ( rows() < nelements ) 
  {
        std::cerr << "CCfits: More data requested than contained in table. ";
        std::cerr << "Extracting complete column.\n";
        nelements = rows();
   }   

   int   status(0);
   int   anynul(0);

   FITSUtil::auto_array_ptr<T> array(new T[nelements]); 

   makeHDUCurrent();

   if ( fits_read_col(fitsPointer(),type(),  index(), firstRow, 1, 
  	nelements, nullValue, array.get(), &anynul, &status) ) throw FitsError(status);


   if (m_data.size() != static_cast<size_t>( rows() ) ) m_data.resize(rows());

   std::copy(&array[0],&array[nelements],m_data.begin()+firstRow-1);
   if (nelements == rows()) isRead(true); 
  }

  template <typename T>
  bool ColumnData<T>::compare (const Column &right) const
  {
  if ( !Column::compare(right) ) return false;
  const ColumnData<T>& that = static_cast<const ColumnData<T>&>(right);
  unsigned int n = m_data.size();
  if ( that.m_data.size() != n ) return false;
  for (unsigned int i = 0; i < n ; i++)
  {
        if (m_data[i] != that.m_data[i]) return false;   
  }
  return true;
  }

  template <typename T>
  ColumnData<T>* ColumnData<T>::clone () const
  {
        return new ColumnData<T>(*this);
  }

  template <typename T>
  std::ostream& ColumnData<T>::put (std::ostream& s) const
  {
  Column::put(s);
  if (FITS::verboseMode() && type() != Tstring)
  {
        s << " Column Legal limits: ( " << m_minLegalValue << "," << m_maxLegalValue << " )\n" 
        << " Column Data  limits: ( " << m_minDataValue << "," << m_maxDataValue << " )\n";
  }
  if (!m_data.empty())
  {
        std::ostream_iterator<T> output(s,"\n");
        // output each row on a separate line.
        // user can supply manipulators to stream for formatting.
        std::copy(m_data.begin(),m_data.end(),output);
  }

    return s;
  }

  template <typename T>
  void ColumnData<T>::writeData (T* indata, long nRows, long firstRow, T* nullValue)
  {

          // set columnData's data member to equal what's written to file.
          // indata has size nRows: elements firstRow to firstRow + nRows - 1 will be written.
          // if this exceeds the current rowlength of the HDU, update the return value for
          // rows() in the parent after the fitsio call.
          int status(0);
          long elementsToWrite(nRows + firstRow -1);
          // get a copy for restorative action.   
          std::vector<T> __tmp(m_data);


          if (elementsToWrite > static_cast<long>(m_data.size())) 
          {

                  m_data.resize(elementsToWrite,T());
          }

          std::copy(&indata[0],&indata[nRows],m_data.begin()+firstRow-1);

          // if successful, write to disk.

          try
          {
             if (nullValue)
             {
                if (fits_write_colnull(fitsPointer(), type(), index(), firstRow, 1, nRows,
			          indata, nullValue, &status) != 0) throw FitsError(status);
             }
             else
             {
                if (fits_write_col(fitsPointer(), type(), index(), firstRow, 1, nRows,
			          indata, &status) != 0) throw FitsError(status);
             }

                // tell the Table that the number of rows has changed
                parent()->updateRows();
          }
          catch (FitsError) // the only thing that can throw here.
          {
                  // reset to original content and rethrow the exception.
                  m_data = __tmp;
                  if (status == NO_NULL) throw NoNullValue(name());
                  else throw;
          }      
  }

  template <typename T>
  void ColumnData<T>::writeData (const std::vector<T>& indata, long firstRow, T* nullValue)
  {
        FITSUtil::CVarray<T> convert;
        FITSUtil::auto_array_ptr<T> pcolData (convert(indata));
        T* columnData  = pcolData.get();
        writeData(columnData,indata.size(),firstRow,nullValue);
  }

  template <typename T>
  void ColumnData<T>::insertRows (long first, long number)
  {
    FITSUtil::FitsNullValue<T> blank;
    typename std::vector<T>::iterator in;
    if (first !=0) 
    {
            in = m_data.begin()+first;
    }
    else
    {
            in = m_data.begin();
    }           

    // non-throwing operations.
    m_data.insert(in,number,blank());
  }

  template <typename T>
  void ColumnData<T>::deleteRows (long first, long number)
  {
    m_data.erase(m_data.begin()+first-1,m_data.begin()+first-1+number);
  }

  template <typename T>
  void ColumnData<T>::setDataLimits (T* limits)
  {
    m_minLegalValue = limits[0];
    m_maxLegalValue = limits[1];
    m_minDataValue = std::max(limits[2],limits[0]);
    m_maxDataValue = std::min(limits[3],limits[1]);
  }

  // Additional Declarations

  // all functions that operate on strings or complex data that call cfitsio 
  // need to be specialized.

#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
template <>
inline void ColumnData<complex<float> >::setDataLimits (complex<float>* limits)
       {
                m_minLegalValue = limits[0];
                m_maxLegalValue = limits[1];
                m_minDataValue =  limits[2];
                m_maxDataValue =  limits[3];
        }
#else
template <>
  void ColumnData<complex<float> >::setDataLimits (complex<float>* limits);
#endif

#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
template <>
inline void ColumnData<complex<double> >::setDataLimits (complex<double>* limits)
        {
                m_minLegalValue = limits[0];
                m_maxLegalValue = limits[1];
                m_minDataValue =  limits[2];
                m_maxDataValue =  limits[3];
        }
#else
 template <>
  void ColumnData<complex<double> >::setDataLimits (complex<double>* limits);
#endif


#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
        template <>
        inline void ColumnData<string>::readColumnData (long firstRow, 
                                        long nelements, 
                                        string* nullValue)
        {
           // nelements = nrows to read.
           if (nelements < 1)
              throw Column::InvalidNumberOfRows((int)nelements);
           if (firstRow < 1 || (firstRow+nelements-1)>rows())
              throw Column::InvalidRowNumber(name());
           
           int status = 0;
           int   anynul = 0;
           
           char** array = new char*[nelements];
           // Initialize pointers to NULL so we can safely delete
           //  during error handling, even if they haven't been allocated.           
           for (long i=0; i<nelements; ++i)
              array[i]=static_cast<char*>(0);
           bool isError = false;

           // Strings are unusual.  The variable length case is still
           //  handled by a ColumnData class, not a ColumnVectorData.
           char* nulval = 0;
           if (nullValue) 
           {
              nulval = const_cast<char*>(nullValue->c_str());
           }
           else
           {
              nulval = new char;
              *nulval = '\0';       
           }
           makeHDUCurrent();
           if (varLength())
           {
              long* strLengths = new long[nelements];
              long* offsets = new long[nelements];
              if (fits_read_descripts(fitsPointer(), index(), firstRow,
                   nelements, strLengths, offsets, &status))
              {
                 isError = true;
              }
              else
              {
                 // For variable length cols, must read 1 and only 1 row
                 //  at a time into array.
                 for (long j=0; j<nelements; ++j)
                 {
                    array[j] = new char[strLengths[j] + 1];
                 }
                 
                 const long lastRow = firstRow+nelements-1;
                 for (long iRow=firstRow; !isError && iRow<=lastRow; ++iRow)
                 {
                    if (fits_read_col_str(fitsPointer(),index(), iRow, 1, 1,
                      nulval, &array[iRow-firstRow], &anynul,&status) )
                       isError=true;
                 }
              }
              delete [] strLengths;
              delete [] offsets;              
           }
           else
           {
              // Fixed length strings, length is stored in Column's m_width.
              for (long j=0; j<nelements; ++j)
              {
                  array[j] = new char[width() + 1];
              }
              if (fits_read_col_str(fitsPointer(),index(), firstRow,1,nelements,
                nulval,array, &anynul,&status))
                isError=true;
           }
           
           if (isError)
           {
              // It's OK to do this even if error occurred before
              //  array rows were allocated.  In that case their pointers
              //  were set to NULL.
              for (long j = 0; j < nelements; ++j)
              {
                 delete [] array[j];
              }     
              delete [] array; 
              delete nulval;
              throw FitsError(status); 
           }

          if (m_data.size() != static_cast<size_t>(rows()))
              setData(std::vector<String>(rows(),String(nulval)));

          for (long j=0; j<nelements; ++j)
          {
             m_data[j - 1 + firstRow] = String(array[j]);
          }

          for (long j=0; j<nelements; j++)
          {
             delete [] array[j];
          }     
          delete [] array; 
          delete nulval; 
          if (nelements == rows()) isRead(true); 

        }
#else 
 template <>
void ColumnData<string>::readColumnData (long firstRow, long nelements, string* nullValue);
#endif


#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
        template <>
        inline void ColumnData<complex<float> >::readColumnData (long firstRow,
                                                long nelements,
                                                complex<float>* nullValue)
        {
          // specialization for ColumnData<string> 
          int status(0);
          int   anynul(0);
          FITSUtil::auto_array_ptr<float> pArray(new float[nelements*2]); 
          float* array = pArray.get();
          float nulval(0);
          makeHDUCurrent();


          if (fits_read_col_cmp(fitsPointer(),index(), firstRow,1,nelements,
                  nulval,array, &anynul,&status) ) throw FitsError(status);


          if (m_data.size() != rows()) m_data.resize(rows());

          // the 'j -1 ' converts to zero based indexing.

          for (int j = 0; j < nelements; ++j)
          {

                m_data[j - 1 + firstRow] = std::complex<float>(array[2*j],array[2*j+1]);
          }
          if (nelements == rows()) isRead(true); 

        }
#else
template <> 
void ColumnData<complex<float> >::readColumnData (long firstRow, long nelements,complex<float>* nullValue );
#endif

#if SPEC_TEMPLATE_IMP_DEFECT || SPEC_TEMPLATE_DECL_DEFECT
        template <>
        inline void ColumnData<complex<double> >::readColumnData (long firstRow, 
                                                        long nelements,
                                                        complex<double>* nullValue)
        {
          // specialization for ColumnData<complex<double> > 
           int status(0);
           int   anynul(0);
           FITSUtil::auto_array_ptr<double> pArray(new double[nelements*2]); 
           double* array = pArray.get();
           double nulval(0);
           makeHDUCurrent();


          if (fits_read_col_dblcmp(fitsPointer(), index(), firstRow,1,nelements,
                  nulval,array, &anynul,&status) ) throw FitsError(status);




          if (m_data.size() != rows()) setData(std::vector<complex<double> >(rows(),nulval));

          // the 'j -1 ' converts to zero based indexing.

          for (int j = 0; j < nelements; j++)
          {

                m_data[j - 1 + firstRow] = std::complex<double>(array[2*j],array[2*j+1]);
          }
          if (nelements == rows()) isRead(true); 

        }
#else
template <>
void ColumnData<complex<double> >::readColumnData (long firstRow, long nelements,complex<double>* nullValue);
#endif

#if SPEC_TEMPLATE_DECL_DEFECT
  template <>
  inline void ColumnData<string>::writeData (const std::vector<string>& indata, 
				      long firstRow, string* nullValue)
  {
    int    status=0;
    char** columnData=FITSUtil::CharArray(indata);

    if ( fits_write_colnull(fitsPointer(), TSTRING, index(), firstRow, 1, indata.size(),
			    columnData, 0, &status) != 0 )
      throw FitsError(status);
    unsigned long elementsToWrite (indata.size() + firstRow - 1);
    std::vector<string> __tmp(m_data);
    if (m_data.size() < elementsToWrite) 
      {
	m_data.resize(elementsToWrite,"");
	std::copy(__tmp.begin(),__tmp.end(),m_data.begin());
      }
    std::copy(indata.begin(),indata.end(),m_data.begin()+firstRow-1);


    for (size_t i = 0; i < indata.size(); ++i)
      {
	delete [] columnData[i];
      }
    delete [] columnData;
  }  
#else
template <>
void ColumnData<string>::writeData (const std::vector<string>& inData, long firstRow, string* nullValue);
#endif

#ifdef SPEC_TEMPLATE_DECL_DEFECT
  template <>
  inline void ColumnData<complex<float> >::writeData (const std::vector<complex<float> >& inData, 
					       long firstRow, 
					       complex<float>* nullValue)
  {
    int status(0);
    int nRows (inData.size());
    FITSUtil::auto_array_ptr<float> pData(new float[nRows*2]);
    float* Data = pData.get();
    std::vector<complex<float> > __tmp(m_data);
    for (int j = 0; j < nRows; ++j)
      {
	Data[ 2*j] = inData[j].real();
	Data[ 2*j + 1] = inData[j].imag();
      }     

    try
      {

	if (fits_write_col_cmp(fitsPointer(), index(), firstRow, 1, 
			       nRows,Data, &status) != 0) throw FitsError(status);
	long elementsToWrite(nRows + firstRow -1);
	if (elementsToWrite > static_cast<long>(m_data.size())) 
	  {

	    m_data.resize(elementsToWrite);
	  }

	std::copy(inData.begin(),inData.end(),m_data.begin()+firstRow-1);

	// tell the Table that the number of rows has changed
	parent()->updateRows();
      }
    catch (FitsError) // the only thing that can throw here.
      {
	// reset to original content and rethrow the exception.
	m_data.resize(__tmp.size());
	m_data = __tmp;
      }      

  }

#else
template <>
void ColumnData<complex<float> >::writeData (const std::vector<complex<float> >& inData, long firstRow, 
                                complex<float>* nullValue);
#endif

#ifdef SPEC_TEMPLATE_DECL_DEFECT
  template <>
  inline void ColumnData<complex<double> >::writeData (const std::vector<complex<double> >& inData, 
						long firstRow, 
						complex<double>* nullValue)
  {
    int status(0);
    int nRows (inData.size());
    FITSUtil::auto_array_ptr<double> pData(new double[nRows*2]);
    double* Data = pData.get();
    std::vector<complex<double> > __tmp(m_data);
    for (int j = 0; j < nRows; ++j)
      {
	pData[ 2*j] = inData[j].real();
	pData[ 2*j + 1] = inData[j].imag();
      }     

    try
      {

	if (fits_write_col_dblcmp(fitsPointer(), index(), firstRow, 1, 
                                  nRows,Data, &status) != 0) throw FitsError(status);
	long elementsToWrite(nRows + firstRow -1);
	if (elementsToWrite > static_cast<long>(m_data.size())) 
	  {

	    m_data.resize(elementsToWrite);
	  }

	std::copy(inData.begin(),inData.end(),m_data.begin()+firstRow-1);

	// tell the Table that the number of rows has changed
	parent()->updateRows();
      }
    catch (FitsError) // the only thing that can throw here.
      {
	// reset to original content and rethrow the exception.
	m_data.resize(__tmp.size());
	m_data = __tmp;
      }      

  }

#else
template <>
void ColumnData<complex<double> >::writeData (const std::vector<complex<double> >& inData, long firstRow, 
                                complex<double>* nullValue);

#endif
} // namespace CCfits


#endif