source: trunk/src/Scantable.h @ 2290

Last change on this file since 2290 was 2290, checked in by Kana Sugimoto, 13 years ago

New Development: Yes

JIRA Issue: Yes (CAS-3219/ASAP-247)

Ready for Test: Yes

Interface Changes: No

What Interface Changed:

Test Programs: sdlist unit test (reference data to be updated)

Put in Release Notes: Yes

Module(s): scantable.summary, sdlist

Description:

Output format of scantable summary changed.
Less use of TableIterator? for speed up scantable.summary/sdlist now lists a scantable
with 348,000 records (NRO 45m w/ 25beams x 13,920scans) in ~30sec (was ~7 min previously).


  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 23.9 KB
Line 
1//
2// C++ Interface: Scantable
3//
4// Description:
5//
6//
7// Author: Malte Marquarding <asap@atnf.csiro.au>, (C) 2005
8//
9// Copyright: See COPYING file that comes with this distribution
10//
11//
12#ifndef ASAPSCANTABLE_H
13#define ASAPSCANTABLE_H
14
15// STL
16#include <fstream>
17#include <iostream>
18#include <sstream>
19#include <string>
20#include <vector>
21// AIPS++
22#include <casa/aips.h>
23#include <casa/Arrays/MaskedArray.h>
24#include <casa/Arrays/Vector.h>
25#include <casa/BasicSL/String.h>
26#include <casa/Containers/Record.h>
27#include <casa/Exceptions/Error.h>
28#include <casa/Quanta/Quantum.h>
29#include <casa/Utilities/CountedPtr.h>
30
31#include <coordinates/Coordinates/SpectralCoordinate.h>
32
33#include <measures/TableMeasures/ScalarMeasColumn.h>
34
35#include <scimath/Mathematics/FFTServer.h>
36
37#include <tables/Tables/ArrayColumn.h>
38#include <tables/Tables/ScalarColumn.h>
39#include <tables/Tables/Table.h>
40
41#include "Logger.h"
42#include "MathUtils.h"
43#include "STFit.h"
44#include "STFitEntry.h"
45#include "STFitter.h"
46#include "STFocus.h"
47#include "STFrequencies.h"
48#include "STHeader.h"
49#include "STHistory.h"
50#include "STMolecules.h"
51#include "STPol.h"
52#include "STSelector.h"
53#include "STTcal.h"
54#include "STWeather.h"
55
56namespace asap {
57
58/**
59  * This class contains and wraps a casa::Table, which is used to store
60  * all the information. This can be either a MemoryTable or a
61  * disk based Table.
62  * It provides access functions to the underlying table
63  * It contains n subtables:
64  * @li weather
65  * @li frequencies
66  * @li molecules
67  * @li tcal
68  * @li focus
69  * @li fits
70  *
71  * @brief The main ASAP data container
72  * @author Malte Marquarding
73  * @date
74  * @version
75*/
76class Scantable : private Logger
77{
78
79friend class STMath;
80
81public:
82  /**
83   * Default constructor
84   */
85  explicit Scantable(casa::Table::TableType ttype = casa::Table::Memory);
86
87  /**
88   * Create a Scantable object form an existing table on disk
89   * @param[in] name the name of the existing Scantable
90   */
91  explicit Scantable(const std::string& name, casa::Table::TableType ttype = casa::Table::Memory);
92
93  /// @fixme this is only sensible for MemoryTables....
94  Scantable(const Scantable& other, bool clear=true);
95
96  /**
97   * Destructor
98   */
99  virtual ~Scantable();
100
101  /**
102   * get a const reference to the underlying casa::Table
103   * @return const \ref casa::Table reference
104   */
105  const casa::Table& table() const;
106
107  /**
108   * get a reference to the underlying casa::Table with the Selection
109   * object applied if set
110   * @return casa::Table reference
111   */
112  casa::Table& table();
113
114
115  /**
116   * Get a handle to the selection object
117   * @return constant STSelector reference
118   */
119  const STSelector& getSelection() const { return selector_; }
120
121  /**
122   * Set the data to be a subset as defined by the STSelector
123   * @param selection a STSelector object
124   */
125  void setSelection(const STSelector& selection);
126
127  /**
128   * unset the selection of the data
129   */
130  void unsetSelection();
131  /**
132   * set the header
133   * @param[in] sth an STHeader object
134   */
135  void setHeader( const STHeader& sth );
136
137  /**
138   * get the header information
139   * @return an STHeader object
140   */
141  STHeader getHeader( ) const;
142
143  /**
144   * Checks if the "other" Scantable is conformant with this,
145   * i.e. if  header values are the same.
146   * @param[in] other another Scantable
147   * @return true or false
148   */
149  bool conformant( const Scantable& other);
150
151  /**
152   *
153   * @param stype The type of the source, 0 = on, 1 = off
154   */
155  void setSourceType(int stype);
156
157
158  /**
159   * The number of scans in the table
160   * @return number of scans in the table
161   */
162  int nscan() const;
163
164  casa::MEpoch::Types getTimeReference() const;
165
166
167  casa::MEpoch getEpoch(int whichrow) const;
168
169  /**
170   * Get global antenna position
171   * @return casa::MPosition
172   */
173  casa::MPosition getAntennaPosition() const;
174
175  /**
176   * the @ref casa::MDirection for a specific row
177   * @param[in] whichrow the row number
178   * return casa::MDirection
179   */
180  casa::MDirection getDirection( int whichrow ) const;
181
182  /**
183   * get the direction type as a string, e.g. "J2000"
184   * @param[in] whichrow the row number
185   * return the direction string
186   */
187  std::string getDirectionString( int whichrow ) const;
188
189  /**
190   * set the direction type as a string, e.g. "J2000"
191   * @param[in] refstr the direction type
192   */
193  void setDirectionRefString(const std::string& refstr="");
194
195  /**
196   * get the direction reference string
197   * @return a string describing the direction reference
198   */
199  std::string getDirectionRefString() const;
200
201  /**
202   *  Return the Flux unit of the data, e.g. "Jy" or "K"
203   * @return string
204   */
205  std::string getFluxUnit() const;
206
207  /**
208   * Set the Flux unit of the data
209   * @param unit a string representing the unit, e.g "Jy" or "K"
210   */
211  void setFluxUnit( const std::string& unit );
212
213  /**
214   * Set the Stokes type of the data
215   * @param feedtype a string representing the type, e.g "circular" or "linear"
216   */
217  void setFeedType( const std::string& feedtype );
218
219  /**
220   *
221   * @param instrument a string representing an insturment. see xxx
222   */
223  void setInstrument( const std::string& instrument );
224
225  /**
226   * (Re)calculate the azimuth and elevationnfor all rows
227   */
228  void calculateAZEL();
229
230  /**
231   * "hard" flag the data, this flags everything selected in setSelection()
232   * param[in] msk a boolean mask of length nchan describing the points to
233   * to be flagged
234   */
235  //void flag( const std::vector<bool>& msk = std::vector<bool>());
236  //void flag( const std::vector<bool>& msk = std::vector<bool>(), bool unflag=false);
237
238  void flag( int whichrow = -1, const std::vector<bool>& msk = std::vector<bool>(), bool unflag=false);
239
240  /**
241   * Flag the data in a row-based manner. (CAS-1433 Wataru Kawasaki)
242   * param[in] rows    list of row numbers to be flagged
243   */
244  void flagRow( const std::vector<casa::uInt>& rows = std::vector<casa::uInt>(), bool unflag=false);
245
246  /**
247   * Get flagRow info at the specified row. If true, the whole data
248   * at the row should be flagged.
249   */
250  bool getFlagRow(int whichrow) const
251    { return (flagrowCol_(whichrow) > 0); }
252
253  /**
254   * Flag the data outside a specified range (in a channel-based manner).
255   * (CAS-1807 Wataru Kawasaki)
256   */
257  void clip(const casa::Float uthres, const casa::Float dthres, bool clipoutside, bool unflag);
258
259  /**
260   * Return a list of booleans with the size of nchan for a specified row, to get info
261   * about which channel is clipped.
262   */
263  std::vector<bool> getClipMask(int whichrow, const casa::Float uthres, const casa::Float dthres, bool clipoutside, bool unflag);
264  void srchChannelsToClip(casa::uInt whichrow, const casa::Float uthres, const casa::Float dthres, bool clipoutside, bool unflag,
265                          casa::Vector<casa::uChar> flgs);
266
267  /**
268   * Return a list of row numbers with respect to the original table.
269   * @return a list of unsigned ints
270   */
271  std::vector<unsigned int> rownumbers() const;
272
273
274  /**
275   * Get the number of beams in the data or a specific scan
276   * @param scanno the scan number to get the number of beams for.
277   * If scanno<0 the number is retrieved from the header.
278   * @return an integer number
279   */
280  int nbeam(int scanno=-1) const;
281  /**
282   * Get the number of IFs in the data or a specific scan
283   * @param scanno the scan number to get the number of IFs for.
284   * If scanno<0 the number is retrieved from the header.
285   * @return an integer number
286   */
287  int nif(int scanno=-1) const;
288  /**
289   * Get the number of polarizations in the data or a specific scan
290   * @param scanno the scan number to get the number of polarizations for.
291   * If scanno<0 the number is retrieved from the header.
292   * @return an integer number
293   */
294  int npol(int scanno=-1) const;
295
296  std::string getPolType() const;
297
298
299  /**
300   * Get the number of integartion cycles
301   * @param scanno the scan number to get the number of rows for.
302   * If scanno<0 the number is retrieved from the header.
303   * @return the number of rows (for the specified scanno)
304   */
305  int nrow(int scanno=-1) const;
306
307  int getBeam(int whichrow) const;
308  std::vector<uint> getBeamNos() const { return getNumbers(beamCol_); }
309
310  int getIF(int whichrow) const;
311  std::vector<uint> getIFNos() const { return getNumbers(ifCol_); }
312
313  int getPol(int whichrow) const;
314  std::vector<uint> getPolNos() const { return getNumbers(polCol_); }
315
316  std::vector<uint> getScanNos() const { return getNumbers(scanCol_); }
317  int getScan(int whichrow) const { return scanCol_(whichrow); }
318
319  //TT addition
320  std::vector<uint> getMolNos() {return getNumbers(mmolidCol_); }
321
322  /**
323   * Get the number of channels in the data or a specific IF. This currently
324   * varies only with IF number
325   * @param ifno the IF number to get the number of channels for.
326   * If ifno<0 the number is retrieved from the header.
327   * @return an integer number
328   */
329  int nchan(int ifno=-1) const;
330  int getChannels(int whichrow) const;
331
332  int ncycle(int scanno=-1) const;
333  int getCycle(int whichrow) const { return cycleCol_(whichrow); }
334
335  double getInterval(int whichrow) const
336    { return integrCol_(whichrow); }
337
338  float getTsys(int whichrow) const
339    { return casa::Vector<casa::Float>(tsysCol_(whichrow))(0); }
340  std::vector<float> getTsysSpectrum(int whichrow) const ;
341  float getElevation(int whichrow) const
342    { return elCol_(whichrow); }
343  float getAzimuth(int whichrow) const
344    { return azCol_(whichrow); }
345  float getParAngle(int whichrow) const
346    { return focus().getParAngle(mfocusidCol_(whichrow)); }
347  int getTcalId(int whichrow) const
348    { return mtcalidCol_(whichrow); }
349
350  std::string getSourceName(int whichrow) const
351    { return srcnCol_(whichrow); }
352
353  std::vector<bool> getMask(int whichrow) const;
354  std::vector<float> getSpectrum(int whichrow,
355                                 const std::string& poltype = "" ) const;
356
357  void setSpectrum(const std::vector<float>& spec, int whichrow);
358
359  std::string getPolarizationLabel(int index, const std::string& ptype) const
360    { return STPol::getPolLabel(index, ptype ); }
361
362  /**
363   * Write the Scantable to disk
364   * @param filename the output file name
365   */
366  void makePersistent(const std::string& filename);
367
368  std::vector<std::string> getHistory() const
369    { return historyTable_.getHistory(); };
370
371  void addHistory(const std::string& hist) { historyTable_.addEntry(hist); }
372
373  void appendToHistoryTable(const STHistory& otherhist)
374    { historyTable_.append(otherhist); }
375
376  std::string headerSummary();
377  void summary(const std::string& filename="");
378  std::string oldheaderSummary();
379  //  std::string summary();
380  void oldsummary(const std::string& filename="");
381
382  //std::string getTime(int whichrow=-1, bool showdate=true) const;
383  std::string getTime(int whichrow=-1, bool showdate=true, casa::uInt prec=0) const;
384  double getIntTime(int whichrow) const { return integrCol_(whichrow); }
385
386  // returns unit, conversion frame, doppler, base-frame
387
388  /**
389   * Get the frequency set up
390   * This is forwarded to the STFrequencies subtable
391   * @return unit, frame, doppler
392   */
393  std::vector<std::string> getCoordInfo() const
394    { return freqTable_.getInfo(); };
395
396  void setCoordInfo(std::vector<string> theinfo)
397    { return freqTable_.setInfo(theinfo); };
398
399
400  std::vector<double> getAbcissa(int whichrow) const;
401
402  std::vector<float> getWeather(int whichrow) const;
403
404  std::string getAbcissaLabel(int whichrow) const;
405  std::vector<double> getRestFrequencies() const
406    { return moleculeTable_.getRestFrequencies(); }
407  std::vector<double> getRestFrequency(int id) const
408    { return moleculeTable_.getRestFrequency(id); }
409
410  /**
411  void setRestFrequencies(double rf, const std::string& name = "",
412                          const std::string& = "Hz");
413  **/
414  // Modified by Takeshi Nakazato 05/09/2008
415  /***
416  void setRestFrequencies(vector<double> rf, const vector<std::string>& name = "",
417                          const std::string& = "Hz");
418  ***/
419  void setRestFrequencies(vector<double> rf,
420                          const vector<std::string>& name = vector<std::string>(1,""),
421                          const std::string& = "Hz");
422
423  //void setRestFrequencies(const std::string& name);
424  void setRestFrequencies(const vector<std::string>& name);
425
426  void shift(int npix);
427
428  casa::SpectralCoordinate getSpectralCoordinate(int whichrow) const;
429
430  void convertDirection(const std::string& newframe);
431
432  STFrequencies& frequencies() { return freqTable_; }
433  const STFrequencies& frequencies() const { return freqTable_; }
434  STWeather& weather() { return weatherTable_; }
435  const STWeather& weather() const { return weatherTable_; }
436  STFocus& focus() { return focusTable_; }
437  const STFocus& focus() const { return focusTable_; }
438  STTcal& tcal() { return tcalTable_; }
439  const STTcal& tcal() const { return tcalTable_; }
440  STMolecules& molecules() { return moleculeTable_; }
441  const STMolecules& molecules() const { return moleculeTable_; }
442  STHistory& history() { return historyTable_; }
443  const STHistory& history() const { return historyTable_; }
444  STFit& fit() { return fitTable_; }
445  const STFit& fit() const { return fitTable_; }
446
447  std::vector<std::string> columnNames() const;
448
449  void addFit(const STFitEntry& fit, int row);
450  STFitEntry getFit(int row) const
451    { STFitEntry fe; fitTable_.getEntry(fe, mfitidCol_(row)); return fe; }
452
453  //Added by TT
454  /**
455   * Get the antenna name
456   * @return antenna name string
457   */
458  casa::String getAntennaName() const;
459
460  /**
461   * For GBT MS data only. check a scan list
462   * against the information found in GBT_GO table for
463   * scan number orders to get correct pairs.
464   * @param[in] scan list
465   * @return status
466   */
467  int checkScanInfo(const std::vector<int>& scanlist) const;
468
469  /**
470   * Get the direction as a vector, for a specific row
471   * @param[in] whichrow the row numbyyer
472   * @return the direction in a vector
473   */
474  std::vector<double> getDirectionVector(int whichrow) const;
475
476  /**
477   * Set a flag indicating whether the data was parallactified
478   * @param[in] flag true or false
479   */
480  void parallactify(bool flag)
481    { focus().setParallactify(flag); }
482
483  /**
484   * Reshape spectrum
485   * @param[in] nmin, nmax minimum and maximum channel
486   * @param[in] irow       row number
487   *
488   * 30/07/2008 Takeshi Nakazato
489   **/
490  void reshapeSpectrum( int nmin, int nmax ) throw( casa::AipsError );
491  void reshapeSpectrum( int nmin, int nmax, int irow ) ;
492
493  /**
494   * Change channel number under fixed bandwidth
495   * @param[in] nchan, dnu new channel number and spectral resolution
496   * @param[in] irow       row number
497   *
498   * 27/08/2008 Takeshi Nakazato
499   **/
500  void regridChannel( int nchan, double dnu ) ;
501  void regridChannel( int nchan, double dnu, int irow ) ;
502
503  bool getFlagtraFast(casa::uInt whichrow);
504
505  void polyBaseline(const std::vector<bool>& mask,
506                    int order,
507                    bool getResidual=true,
508                    const std::string& progressInfo="true,1000",
509                    const bool outLogger=false,
510                    const std::string& blfile="");
511  void autoPolyBaseline(const std::vector<bool>& mask,
512                        int order,
513                        const std::vector<int>& edge,
514                        float threshold=3.0,
515                        int chanAvgLimit=1,
516                        bool getResidual=true,
517                        const std::string& progressInfo="true,1000",
518                        const bool outLogger=false,
519                        const std::string& blfile="");
520  void cubicSplineBaseline(const std::vector<bool>& mask,
521                           int nPiece,
522                           float thresClip,
523                           int nIterClip,
524                           bool getResidual=true,
525                           const std::string& progressInfo="true,1000",
526                           const bool outLogger=false,
527                           const std::string& blfile="");
528  void autoCubicSplineBaseline(const std::vector<bool>& mask,
529                               int nPiece,
530                               float thresClip,
531                               int nIterClip,
532                               const std::vector<int>& edge,
533                               float threshold=3.0,
534                               int chanAvgLimit=1,
535                               bool getResidual=true,
536                               const std::string& progressInfo="true,1000",
537                               const bool outLogger=false,
538                               const std::string& blfile="");
539  void sinusoidBaseline(const std::vector<bool>& mask,
540                        const bool applyFFT,
541                        const std::string& fftMethod,
542                        const std::string& fftThresh,
543                        const std::vector<int>& addNWaves,
544                        const std::vector<int>& rejectNWaves,
545                        float thresClip,
546                        int nIterClip,
547                        bool getResidual=true,
548                        const std::string& progressInfo="true,1000",
549                        const bool outLogger=false,
550                        const std::string& blfile="");
551  void autoSinusoidBaseline(const std::vector<bool>& mask,
552                            const bool applyFFT,
553                            const std::string& fftMethod,
554                            const std::string& fftThresh,
555                            const std::vector<int>& addNWaves,
556                            const std::vector<int>& rejectNWaves,
557                            float thresClip,
558                            int nIterClip,
559                            const std::vector<int>& edge,
560                            float threshold=3.0,
561                            int chanAvgLimit=1,
562                            bool getResidual=true,
563                            const std::string& progressInfo="true,1000",
564                            const bool outLogger=false,
565                            const std::string& blfile="");
566  std::vector<float> execFFT(const int whichrow,
567                             const std::vector<bool>& inMask,
568                             bool getRealImag=false,
569                             bool getAmplitudeOnly=false);
570  float getRms(const std::vector<bool>& mask, int whichrow);
571  std::string formatBaselineParams(const std::vector<float>& params,
572                                   const std::vector<bool>& fixed,
573                                   float rms,
574                                   int nClipped,
575                                   const std::string& masklist,
576                                   int whichrow,
577                                   bool verbose=false,
578                                   int start=-1,
579                                   int count=-1,
580                                   bool resetparamid=false) const;
581  std::string formatPiecewiseBaselineParams(const std::vector<int>& ranges,
582                                            const std::vector<float>& params,
583                                            const std::vector<bool>& fixed,
584                                            float rms,
585                                            int nClipped,
586                                            const std::string& masklist,
587                                            int whichrow,
588                                            bool verbose=false) const;
589
590
591private:
592
593  casa::Matrix<casa::Float> getPolMatrix( casa::uInt whichrow ) const;
594
595  /**
596   * Turns a time value into a formatted string
597   * @param x
598   * @return
599   */
600  std::string formatSec(casa::Double x) const;
601
602  std::string formatTime(const casa::MEpoch& me, bool showdate)const;
603  std::string formatTime(const casa::MEpoch& me, bool showdate, casa::uInt prec)const;
604
605  /**
606   *  Turns a casa::MDirection into a nicely formatted string
607   * @param md an casa::MDirection
608   * @return
609   */
610  std::string formatDirection(const casa::MDirection& md) const;
611
612  /**
613   * Create a unique file name for the paged (temporary) table
614   * @return just the name
615   */
616  static casa::String generateName();
617
618  /**
619   * attach to cached columns
620   */
621  void attach();
622
623  /**
624   * Set up the main casa::Table
625   */
626  void setupMainTable();
627
628  void attachSubtables();
629  void copySubtables(const Scantable& other);
630
631  /**
632   * Convert an "old" asap1 style row index into a new index
633   * @param[in] therow
634   * @return and index into @table_
635   */
636  int rowToScanIndex(int therow);
637
638  std::vector<uint> getNumbers(const casa::ScalarColumn<casa::uInt>& col) const;
639
640  static const casa::uInt version_ = 3;
641
642  STSelector selector_;
643
644  casa::Table::TableType type_;
645
646  // the actual data
647  casa::Table table_;
648  casa::Table originalTable_;
649
650  STTcal tcalTable_;
651  STFrequencies freqTable_;
652  STWeather weatherTable_;
653  STFocus focusTable_;
654  STMolecules moleculeTable_;
655  STHistory historyTable_;
656  STFit fitTable_;
657
658  // Cached Columns to avoid reconstructing them for each row get/put
659  casa::ScalarColumn<casa::Double> integrCol_;
660  casa::MDirection::ScalarColumn dirCol_;
661  casa::MEpoch::ScalarColumn timeCol_;
662  casa::ScalarColumn<casa::Float> azCol_;
663  casa::ScalarColumn<casa::Float> elCol_;
664  casa::ScalarColumn<casa::String> srcnCol_, fldnCol_;
665  casa::ScalarColumn<casa::uInt> scanCol_, beamCol_, ifCol_, polCol_, cycleCol_, flagrowCol_;
666  casa::ScalarColumn<casa::Int> rbeamCol_, srctCol_;
667  casa::ArrayColumn<casa::Float> specCol_, tsysCol_;
668  casa::ArrayColumn<casa::uChar> flagsCol_;
669
670  // id in frequencies table
671  casa::ScalarColumn<casa::uInt> mfreqidCol_;
672  // id in tcal table
673  casa::ScalarColumn<casa::uInt> mtcalidCol_;
674
675  casa::ArrayColumn<casa::String> histitemCol_;
676  casa::ScalarColumn<casa::Int> mfitidCol_;
677  casa::ScalarColumn<casa::uInt> mweatheridCol_;
678
679  casa::ScalarColumn<casa::uInt> mfocusidCol_;
680
681  casa::ScalarColumn<casa::uInt> mmolidCol_;
682
683  static std::map<std::string, STPol::STPolFactory *> factories_;
684  void initFactories();
685
686  /**
687   * Add an auxiliary column to the main table and attach it to a
688   * cached column. Use for adding new columns that the original asap2
689   * tables do not have.
690   * @param[in] col      reference to the cached column to be attached
691   * @param[in] colName  column name in asap table
692   * @param[in] defValue default value to fill in the column
693   *
694   * 25/10/2009 Wataru Kawasaki
695   */
696  template<class T, class T2> void attachAuxColumnDef(casa::ScalarColumn<T>&,
697                                                       const casa::String&,
698                                                       const T2&);
699  template<class T, class T2> void attachAuxColumnDef(casa::ArrayColumn<T>&,
700                                                      const casa::String&,
701                                                      const casa::Array<T2>&);
702
703  void fitBaseline(const std::vector<bool>& mask, int whichrow, Fitter& fitter);
704  std::vector<float> doCubicSplineFitting(const std::vector<float>& data,
705                                          const std::vector<bool>& mask,
706                                          int nPiece,
707                                          std::vector<int>& idxEdge,
708                                          std::vector<float>& params,
709                                          int& nClipped,
710                                          float thresClip=3.0,
711                                          int nIterClip=1,
712                                          bool getResidual=true);
713  std::vector<float> doSinusoidFitting(const std::vector<float>& data,
714                                       const std::vector<bool>& mask,
715                                       const std::vector<int>& waveNumbers,
716                                       std::vector<float>& params,
717                                       int& nClipped,
718                                       float thresClip=3.0,
719                                       int nIterClip=1,
720                                       bool getResidual=true);
721  void selectWaveNumbers(const int whichrow,
722                         const std::vector<bool>& chanMask,
723                         const bool applyFFT,
724                         const std::string& fftMethod,
725                         const std::string& fftThresh,
726                         const std::vector<int>& addNWaves,
727                         const std::vector<int>& rejectNWaves,
728                         std::vector<int>& nWaves);
729  void parseThresholdExpression(const std::string& fftThresh,
730                                std::string& fftThAttr,
731                                float& fftThSigma,
732                                int& fftThTop);
733  void doSelectWaveNumbers(const int whichrow,
734                           const std::vector<bool>& chanMask,
735                           const std::string& fftMethod,
736                           const float fftThSigma,
737                           const int fftThTop,
738                           const std::string& fftThAttr,
739                           std::vector<int>& nWaves);
740  void addAuxWaveNumbers(const std::vector<int>& addNWaves,
741                         const std::vector<int>& rejectNWaves,
742                         std::vector<int>& nWaves);
743  bool hasSameNchanOverIFs();
744  std::string getMaskRangeList(const std::vector<bool>& mask,
745                                int whichrow,
746                                const casa::String& coordInfo,
747                                bool hasSameNchan,
748                                bool verbose=false);
749  std::vector<int> getMaskEdgeIndices(const std::vector<bool>& mask);
750  std::string formatBaselineParamsHeader(int whichrow, const std::string& masklist, bool verbose) const;
751  std::string formatBaselineParamsFooter(float rms, int nClipped, bool verbose) const;
752  std::vector<bool> getCompositeChanMask(int whichrow, const std::vector<bool>& inMask);
753  //std::vector<bool> getCompositeChanMask(int whichrow, const std::vector<bool>& inMask, const std::vector<int>& edge, const int minEdgeSize, STLineFinder& lineFinder);
754  void outputFittingResult(bool outLogger, bool outTextFile, const std::vector<bool>& chanMask, int whichrow, const casa::String& coordInfo, bool hasSameNchan, std::ofstream& ofs, const casa::String& funcName, Fitter& fitter);
755  void outputFittingResult(bool outLogger, bool outTextFile, const std::vector<bool>& chanMask, int whichrow, const casa::String& coordInfo, bool hasSameNchan, std::ofstream& ofs, const casa::String& funcName, const std::vector<int>& edge, const std::vector<float>& params, const int nClipped);
756  void outputFittingResult(bool outLogger, bool outTextFile, const std::vector<bool>& chanMask, int whichrow, const casa::String& coordInfo, bool hasSameNchan, std::ofstream& ofs, const casa::String& funcName, const std::vector<float>& params, const int nClipped);
757  void parseProgressInfo(const std::string& progressInfo, bool& showProgress, int& minNRow);
758  void showProgressOnTerminal(const int nProcessed, const int nTotal, const bool showProgress=true, const int nTotalThreshold=1000);
759
760  void applyChanFlag( casa::uInt whichrow, const std::vector<bool>& msk, casa::uChar flagval);
761
762
763};
764
765} // namespace
766
767#endif
Note: See TracBrowser for help on using the repository browser.