// // C++ Implementation: Scantable // // Description: // // // Author: Malte Marquarding , (C) 2005 // // Copyright: See COPYING file that comes with this distribution // // #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include // needed to avoid error in .tcc #include // #include #include #include #include #include #include #include #include #include #include "Scantable.h" #include "STPolLinear.h" #include "STPolCircular.h" #include "STPolStokes.h" #include "STAttr.h" #include "STLineFinder.h" #include "MathUtils.h" using namespace casa; namespace asap { std::map Scantable::factories_; void Scantable::initFactories() { if ( factories_.empty() ) { Scantable::factories_["linear"] = &STPolLinear::myFactory; Scantable::factories_["circular"] = &STPolCircular::myFactory; Scantable::factories_["stokes"] = &STPolStokes::myFactory; } } Scantable::Scantable(Table::TableType ttype) : type_(ttype) { initFactories(); setupMainTable(); freqTable_ = STFrequencies(*this); table_.rwKeywordSet().defineTable("FREQUENCIES", freqTable_.table()); weatherTable_ = STWeather(*this); table_.rwKeywordSet().defineTable("WEATHER", weatherTable_.table()); focusTable_ = STFocus(*this); table_.rwKeywordSet().defineTable("FOCUS", focusTable_.table()); tcalTable_ = STTcal(*this); table_.rwKeywordSet().defineTable("TCAL", tcalTable_.table()); moleculeTable_ = STMolecules(*this); table_.rwKeywordSet().defineTable("MOLECULES", moleculeTable_.table()); historyTable_ = STHistory(*this); table_.rwKeywordSet().defineTable("HISTORY", historyTable_.table()); fitTable_ = STFit(*this); table_.rwKeywordSet().defineTable("FIT", fitTable_.table()); table_.tableInfo().setType( "Scantable" ) ; originalTable_ = table_; attach(); } Scantable::Scantable(const std::string& name, Table::TableType ttype) : type_(ttype) { initFactories(); Table tab(name, Table::Update); uInt version = tab.keywordSet().asuInt("VERSION"); if (version != version_) { throw(AipsError("Unsupported version of ASAP file.")); } if ( type_ == Table::Memory ) { table_ = tab.copyToMemoryTable(generateName()); } else { table_ = tab; } table_.tableInfo().setType( "Scantable" ) ; attachSubtables(); originalTable_ = table_; attach(); } /* Scantable::Scantable(const std::string& name, Table::TableType ttype) : type_(ttype) { initFactories(); Table tab(name, Table::Update); uInt version = tab.keywordSet().asuInt("VERSION"); if (version != version_) { throw(AipsError("Unsupported version of ASAP file.")); } if ( type_ == Table::Memory ) { table_ = tab.copyToMemoryTable(generateName()); } else { table_ = tab; } attachSubtables(); originalTable_ = table_; attach(); } */ Scantable::Scantable( const Scantable& other, bool clear ) { // with or without data String newname = String(generateName()); type_ = other.table_.tableType(); if ( other.table_.tableType() == Table::Memory ) { if ( clear ) { table_ = TableCopy::makeEmptyMemoryTable(newname, other.table_, True); } else table_ = other.table_.copyToMemoryTable(newname); } else { other.table_.deepCopy(newname, Table::New, False, other.table_.endianFormat(), Bool(clear)); table_ = Table(newname, Table::Update); table_.markForDelete(); } table_.tableInfo().setType( "Scantable" ) ; /// @todo reindex SCANNO, recompute nbeam, nif, npol if ( clear ) copySubtables(other); attachSubtables(); originalTable_ = table_; attach(); } void Scantable::copySubtables(const Scantable& other) { Table t = table_.rwKeywordSet().asTable("FREQUENCIES"); TableCopy::copyRows(t, other.freqTable_.table()); t = table_.rwKeywordSet().asTable("FOCUS"); TableCopy::copyRows(t, other.focusTable_.table()); t = table_.rwKeywordSet().asTable("WEATHER"); TableCopy::copyRows(t, other.weatherTable_.table()); t = table_.rwKeywordSet().asTable("TCAL"); TableCopy::copyRows(t, other.tcalTable_.table()); t = table_.rwKeywordSet().asTable("MOLECULES"); TableCopy::copyRows(t, other.moleculeTable_.table()); t = table_.rwKeywordSet().asTable("HISTORY"); TableCopy::copyRows(t, other.historyTable_.table()); t = table_.rwKeywordSet().asTable("FIT"); TableCopy::copyRows(t, other.fitTable_.table()); } void Scantable::attachSubtables() { freqTable_ = STFrequencies(table_); focusTable_ = STFocus(table_); weatherTable_ = STWeather(table_); tcalTable_ = STTcal(table_); moleculeTable_ = STMolecules(table_); historyTable_ = STHistory(table_); fitTable_ = STFit(table_); } Scantable::~Scantable() { //cout << "~Scantable() " << this << endl; } void Scantable::setupMainTable() { TableDesc td("", "1", TableDesc::Scratch); td.comment() = "An ASAP Scantable"; td.rwKeywordSet().define("VERSION", uInt(version_)); // n Cycles td.addColumn(ScalarColumnDesc("SCANNO")); // new index every nBeam x nIF x nPol td.addColumn(ScalarColumnDesc("CYCLENO")); td.addColumn(ScalarColumnDesc("BEAMNO")); td.addColumn(ScalarColumnDesc("IFNO")); // linear, circular, stokes td.rwKeywordSet().define("POLTYPE", String("linear")); td.addColumn(ScalarColumnDesc("POLNO")); td.addColumn(ScalarColumnDesc("FREQ_ID")); td.addColumn(ScalarColumnDesc("MOLECULE_ID")); ScalarColumnDesc refbeamnoColumn("REFBEAMNO"); refbeamnoColumn.setDefault(Int(-1)); td.addColumn(refbeamnoColumn); ScalarColumnDesc flagrowColumn("FLAGROW"); flagrowColumn.setDefault(uInt(0)); td.addColumn(flagrowColumn); td.addColumn(ScalarColumnDesc("TIME")); TableMeasRefDesc measRef(MEpoch::UTC); // UTC as default TableMeasValueDesc measVal(td, "TIME"); TableMeasDesc mepochCol(measVal, measRef); mepochCol.write(td); td.addColumn(ScalarColumnDesc("INTERVAL")); td.addColumn(ScalarColumnDesc("SRCNAME")); // Type of source (on=0, off=1, other=-1) ScalarColumnDesc stypeColumn("SRCTYPE"); stypeColumn.setDefault(Int(-1)); td.addColumn(stypeColumn); td.addColumn(ScalarColumnDesc("FIELDNAME")); //The actual Data Vectors td.addColumn(ArrayColumnDesc("SPECTRA")); td.addColumn(ArrayColumnDesc("FLAGTRA")); td.addColumn(ArrayColumnDesc("TSYS")); td.addColumn(ArrayColumnDesc("DIRECTION", IPosition(1,2), ColumnDesc::Direct)); TableMeasRefDesc mdirRef(MDirection::J2000); // default TableMeasValueDesc tmvdMDir(td, "DIRECTION"); // the TableMeasDesc gives the column a type TableMeasDesc mdirCol(tmvdMDir, mdirRef); // a uder set table type e.g. GALCTIC, B1950 ... td.rwKeywordSet().define("DIRECTIONREF", String("J2000")); // writing create the measure column mdirCol.write(td); td.addColumn(ScalarColumnDesc("AZIMUTH")); td.addColumn(ScalarColumnDesc("ELEVATION")); td.addColumn(ScalarColumnDesc("OPACITY")); td.addColumn(ScalarColumnDesc("TCAL_ID")); ScalarColumnDesc fitColumn("FIT_ID"); fitColumn.setDefault(Int(-1)); td.addColumn(fitColumn); td.addColumn(ScalarColumnDesc("FOCUS_ID")); td.addColumn(ScalarColumnDesc("WEATHER_ID")); // columns which just get dragged along, as they aren't used in asap td.addColumn(ScalarColumnDesc("SRCVELOCITY")); td.addColumn(ArrayColumnDesc("SRCPROPERMOTION")); td.addColumn(ArrayColumnDesc("SRCDIRECTION")); td.addColumn(ArrayColumnDesc("SCANRATE")); td.rwKeywordSet().define("OBSMODE", String("")); // Now create Table SetUp from the description. SetupNewTable aNewTab(generateName(), td, Table::Scratch); table_ = Table(aNewTab, type_, 0); originalTable_ = table_; } void Scantable::attach() { timeCol_.attach(table_, "TIME"); srcnCol_.attach(table_, "SRCNAME"); srctCol_.attach(table_, "SRCTYPE"); specCol_.attach(table_, "SPECTRA"); flagsCol_.attach(table_, "FLAGTRA"); tsysCol_.attach(table_, "TSYS"); cycleCol_.attach(table_,"CYCLENO"); scanCol_.attach(table_, "SCANNO"); beamCol_.attach(table_, "BEAMNO"); ifCol_.attach(table_, "IFNO"); polCol_.attach(table_, "POLNO"); integrCol_.attach(table_, "INTERVAL"); azCol_.attach(table_, "AZIMUTH"); elCol_.attach(table_, "ELEVATION"); dirCol_.attach(table_, "DIRECTION"); fldnCol_.attach(table_, "FIELDNAME"); rbeamCol_.attach(table_, "REFBEAMNO"); mweatheridCol_.attach(table_,"WEATHER_ID"); mfitidCol_.attach(table_,"FIT_ID"); mfreqidCol_.attach(table_, "FREQ_ID"); mtcalidCol_.attach(table_, "TCAL_ID"); mfocusidCol_.attach(table_, "FOCUS_ID"); mmolidCol_.attach(table_, "MOLECULE_ID"); //Add auxiliary column for row-based flagging (CAS-1433 Wataru Kawasaki) attachAuxColumnDef(flagrowCol_, "FLAGROW", 0); } template void Scantable::attachAuxColumnDef(ScalarColumn& col, const String& colName, const T2& defValue) { try { col.attach(table_, colName); } catch (TableError& err) { String errMesg = err.getMesg(); if (errMesg == "Table column " + colName + " is unknown") { table_.addColumn(ScalarColumnDesc(colName)); col.attach(table_, colName); col.fillColumn(static_cast(defValue)); } else { throw; } } catch (...) { throw; } } template void Scantable::attachAuxColumnDef(ArrayColumn& col, const String& colName, const Array& defValue) { try { col.attach(table_, colName); } catch (TableError& err) { String errMesg = err.getMesg(); if (errMesg == "Table column " + colName + " is unknown") { table_.addColumn(ArrayColumnDesc(colName)); col.attach(table_, colName); int size = 0; ArrayIterator& it = defValue.begin(); while (it != defValue.end()) { ++size; ++it; } IPosition ip(1, size); Array& arr(ip); for (int i = 0; i < size; ++i) arr[i] = static_cast(defValue[i]); col.fillColumn(arr); } else { throw; } } catch (...) { throw; } } void Scantable::setHeader(const STHeader& sdh) { table_.rwKeywordSet().define("nIF", sdh.nif); table_.rwKeywordSet().define("nBeam", sdh.nbeam); table_.rwKeywordSet().define("nPol", sdh.npol); table_.rwKeywordSet().define("nChan", sdh.nchan); table_.rwKeywordSet().define("Observer", sdh.observer); table_.rwKeywordSet().define("Project", sdh.project); table_.rwKeywordSet().define("Obstype", sdh.obstype); table_.rwKeywordSet().define("AntennaName", sdh.antennaname); table_.rwKeywordSet().define("AntennaPosition", sdh.antennaposition); table_.rwKeywordSet().define("Equinox", sdh.equinox); table_.rwKeywordSet().define("FreqRefFrame", sdh.freqref); table_.rwKeywordSet().define("FreqRefVal", sdh.reffreq); table_.rwKeywordSet().define("Bandwidth", sdh.bandwidth); table_.rwKeywordSet().define("UTC", sdh.utc); table_.rwKeywordSet().define("FluxUnit", sdh.fluxunit); table_.rwKeywordSet().define("Epoch", sdh.epoch); table_.rwKeywordSet().define("POLTYPE", sdh.poltype); } STHeader Scantable::getHeader() const { STHeader sdh; table_.keywordSet().get("nBeam",sdh.nbeam); table_.keywordSet().get("nIF",sdh.nif); table_.keywordSet().get("nPol",sdh.npol); table_.keywordSet().get("nChan",sdh.nchan); table_.keywordSet().get("Observer", sdh.observer); table_.keywordSet().get("Project", sdh.project); table_.keywordSet().get("Obstype", sdh.obstype); table_.keywordSet().get("AntennaName", sdh.antennaname); table_.keywordSet().get("AntennaPosition", sdh.antennaposition); table_.keywordSet().get("Equinox", sdh.equinox); table_.keywordSet().get("FreqRefFrame", sdh.freqref); table_.keywordSet().get("FreqRefVal", sdh.reffreq); table_.keywordSet().get("Bandwidth", sdh.bandwidth); table_.keywordSet().get("UTC", sdh.utc); table_.keywordSet().get("FluxUnit", sdh.fluxunit); table_.keywordSet().get("Epoch", sdh.epoch); table_.keywordSet().get("POLTYPE", sdh.poltype); return sdh; } void Scantable::setSourceType( int stype ) { if ( stype < 0 || stype > 1 ) throw(AipsError("Illegal sourcetype.")); TableVector tabvec(table_, "SRCTYPE"); tabvec = Int(stype); } bool Scantable::conformant( const Scantable& other ) { return this->getHeader().conformant(other.getHeader()); } std::string Scantable::formatSec(Double x) const { Double xcop = x; MVTime mvt(xcop/24./3600.); // make days if (x < 59.95) return String(" ") + mvt.string(MVTime::TIME_CLEAN_NO_HM, 7)+"s"; else if (x < 3599.95) return String(" ") + mvt.string(MVTime::TIME_CLEAN_NO_H,7)+" "; else { ostringstream oss; oss << setw(2) << std::right << setprecision(1) << mvt.hour(); oss << ":" << mvt.string(MVTime::TIME_CLEAN_NO_H,7) << " "; return String(oss); } }; std::string Scantable::formatDirection(const MDirection& md) const { Vector t = md.getAngle(Unit(String("rad"))).getValue(); Int prec = 7; MVAngle mvLon(t[0]); String sLon = mvLon.string(MVAngle::TIME,prec); uInt tp = md.getRef().getType(); if (tp == MDirection::GALACTIC || tp == MDirection::SUPERGAL ) { sLon = mvLon(0.0).string(MVAngle::ANGLE_CLEAN,prec); } MVAngle mvLat(t[1]); String sLat = mvLat.string(MVAngle::ANGLE+MVAngle::DIG2,prec); return sLon + String(" ") + sLat; } std::string Scantable::getFluxUnit() const { return table_.keywordSet().asString("FluxUnit"); } void Scantable::setFluxUnit(const std::string& unit) { String tmp(unit); Unit tU(tmp); if (tU==Unit("K") || tU==Unit("Jy")) { table_.rwKeywordSet().define(String("FluxUnit"), tmp); } else { throw AipsError("Illegal unit - must be compatible with Jy or K"); } } void Scantable::setInstrument(const std::string& name) { bool throwIt = true; // create an Instrument to see if this is valid STAttr::convertInstrument(name, throwIt); String nameU(name); nameU.upcase(); table_.rwKeywordSet().define(String("AntennaName"), nameU); } void Scantable::setFeedType(const std::string& feedtype) { if ( Scantable::factories_.find(feedtype) == Scantable::factories_.end() ) { std::string msg = "Illegal feed type "+ feedtype; throw(casa::AipsError(msg)); } table_.rwKeywordSet().define(String("POLTYPE"), feedtype); } MPosition Scantable::getAntennaPosition() const { Vector antpos; table_.keywordSet().get("AntennaPosition", antpos); MVPosition mvpos(antpos(0),antpos(1),antpos(2)); return MPosition(mvpos); } void Scantable::makePersistent(const std::string& filename) { String inname(filename); Path path(inname); /// @todo reindex SCANNO, recompute nbeam, nif, npol inname = path.expandedName(); // 2011/03/04 TN // We can comment out this workaround since the essential bug is // fixed in casacore (r20889 in google code). table_.deepCopy(inname, Table::New); // // WORKAROUND !!! for Table bug // // Remove when fixed in casacore // if ( table_.tableType() == Table::Memory && !selector_.empty() ) { // Table tab = table_.copyToMemoryTable(generateName()); // tab.deepCopy(inname, Table::New); // tab.markForDelete(); // // } else { // table_.deepCopy(inname, Table::New); // } } int Scantable::nbeam( int scanno ) const { if ( scanno < 0 ) { Int n; table_.keywordSet().get("nBeam",n); return int(n); } else { // take the first POLNO,IFNO,CYCLENO as nbeam shouldn't vary with these Table t = table_(table_.col("SCANNO") == scanno); ROTableRow row(t); const TableRecord& rec = row.get(0); Table subt = t( t.col("IFNO") == Int(rec.asuInt("IFNO")) && t.col("POLNO") == Int(rec.asuInt("POLNO")) && t.col("CYCLENO") == Int(rec.asuInt("CYCLENO")) ); ROTableVector v(subt, "BEAMNO"); return int(v.nelements()); } return 0; } int Scantable::nif( int scanno ) const { if ( scanno < 0 ) { Int n; table_.keywordSet().get("nIF",n); return int(n); } else { // take the first POLNO,BEAMNO,CYCLENO as nbeam shouldn't vary with these Table t = table_(table_.col("SCANNO") == scanno); ROTableRow row(t); const TableRecord& rec = row.get(0); Table subt = t( t.col("BEAMNO") == Int(rec.asuInt("BEAMNO")) && t.col("POLNO") == Int(rec.asuInt("POLNO")) && t.col("CYCLENO") == Int(rec.asuInt("CYCLENO")) ); if ( subt.nrow() == 0 ) return 0; ROTableVector v(subt, "IFNO"); return int(v.nelements()); } return 0; } int Scantable::npol( int scanno ) const { if ( scanno < 0 ) { Int n; table_.keywordSet().get("nPol",n); return n; } else { // take the first POLNO,IFNO,CYCLENO as nbeam shouldn't vary with these Table t = table_(table_.col("SCANNO") == scanno); ROTableRow row(t); const TableRecord& rec = row.get(0); Table subt = t( t.col("BEAMNO") == Int(rec.asuInt("BEAMNO")) && t.col("IFNO") == Int(rec.asuInt("IFNO")) && t.col("CYCLENO") == Int(rec.asuInt("CYCLENO")) ); if ( subt.nrow() == 0 ) return 0; ROTableVector v(subt, "POLNO"); return int(v.nelements()); } return 0; } int Scantable::ncycle( int scanno ) const { if ( scanno < 0 ) { Block cols(2); cols[0] = "SCANNO"; cols[1] = "CYCLENO"; TableIterator it(table_, cols); int n = 0; while ( !it.pastEnd() ) { ++n; ++it; } return n; } else { Table t = table_(table_.col("SCANNO") == scanno); ROTableRow row(t); const TableRecord& rec = row.get(0); Table subt = t( t.col("BEAMNO") == Int(rec.asuInt("BEAMNO")) && t.col("POLNO") == Int(rec.asuInt("POLNO")) && t.col("IFNO") == Int(rec.asuInt("IFNO")) ); if ( subt.nrow() == 0 ) return 0; return int(subt.nrow()); } return 0; } int Scantable::nrow( int scanno ) const { return int(table_.nrow()); } int Scantable::nchan( int ifno ) const { if ( ifno < 0 ) { Int n; table_.keywordSet().get("nChan",n); return int(n); } else { // take the first SCANNO,POLNO,BEAMNO,CYCLENO as nbeam shouldn't // vary with these Table t = table_(table_.col("IFNO") == ifno); if ( t.nrow() == 0 ) return 0; ROArrayColumn v(t, "SPECTRA"); return v.shape(0)(0); } return 0; } int Scantable::nscan() const { Vector scannos(scanCol_.getColumn()); uInt nout = genSort( scannos, Sort::Ascending, Sort::QuickSort|Sort::NoDuplicates ); return int(nout); } int Scantable::getChannels(int whichrow) const { return specCol_.shape(whichrow)(0); } int Scantable::getBeam(int whichrow) const { return beamCol_(whichrow); } std::vector Scantable::getNumbers(const ScalarColumn& col) const { Vector nos(col.getColumn()); uInt n = genSort( nos, Sort::Ascending, Sort::QuickSort|Sort::NoDuplicates ); nos.resize(n, True); std::vector stlout; nos.tovector(stlout); return stlout; } int Scantable::getIF(int whichrow) const { return ifCol_(whichrow); } int Scantable::getPol(int whichrow) const { return polCol_(whichrow); } std::string Scantable::formatTime(const MEpoch& me, bool showdate) const { return formatTime(me, showdate, 0); } std::string Scantable::formatTime(const MEpoch& me, bool showdate, uInt prec) const { MVTime mvt(me.getValue()); if (showdate) //mvt.setFormat(MVTime::YMD); mvt.setFormat(MVTime::YMD, prec); else //mvt.setFormat(MVTime::TIME); mvt.setFormat(MVTime::TIME, prec); ostringstream oss; oss << mvt; return String(oss); } void Scantable::calculateAZEL() { MPosition mp = getAntennaPosition(); MEpoch::ROScalarColumn timeCol(table_, "TIME"); ostringstream oss; oss << "Computed azimuth/elevation using " << endl << mp << endl; for (Int i=0; i "; MeasFrame frame(mp, me); Vector azel = MDirection::Convert(md, MDirection::Ref(MDirection::AZEL, frame) )().getAngle("rad").getValue(); azCol_.put(i,Float(azel[0])); elCol_.put(i,Float(azel[1])); oss << "azel: " << azel[0]/C::pi*180.0 << " " << azel[1]/C::pi*180.0 << " (deg)" << endl; } pushLog(String(oss)); } void Scantable::clip(const Float uthres, const Float dthres, bool clipoutside, bool unflag) { for (uInt i=0; i flgs = flagsCol_(i); srchChannelsToClip(i, uthres, dthres, clipoutside, unflag, flgs); flagsCol_.put(i, flgs); } } std::vector Scantable::getClipMask(int whichrow, const Float uthres, const Float dthres, bool clipoutside, bool unflag) { Vector flags; flagsCol_.get(uInt(whichrow), flags); srchChannelsToClip(uInt(whichrow), uthres, dthres, clipoutside, unflag, flags); Vector bflag(flags.shape()); convertArray(bflag, flags); //bflag = !bflag; std::vector mask; bflag.tovector(mask); return mask; } void Scantable::srchChannelsToClip(uInt whichrow, const Float uthres, const Float dthres, bool clipoutside, bool unflag, Vector flgs) { Vector spcs = specCol_(whichrow); uInt nchannel = nchan(); if (spcs.nelements() != nchannel) { throw(AipsError("Data has incorrect number of channels")); } uChar userflag = 1 << 7; if (unflag) { userflag = 0 << 7; } if (clipoutside) { for (uInt j = 0; j < nchannel; ++j) { Float spc = spcs(j); if ((spc >= uthres) || (spc <= dthres)) { flgs(j) = userflag; } } } else { for (uInt j = 0; j < nchannel; ++j) { Float spc = spcs(j); if ((spc < uthres) && (spc > dthres)) { flgs(j) = userflag; } } } } void Scantable::flag( int whichrow, const std::vector& msk, bool unflag ) { std::vector::const_iterator it; uInt ntrue = 0; if (whichrow >= int(table_.nrow()) ) { throw(AipsError("Invalid row number")); } for (it = msk.begin(); it != msk.end(); ++it) { if ( *it ) { ntrue++; } } //if ( selector_.empty() && (msk.size() == 0 || msk.size() == ntrue) ) if ( whichrow == -1 && !unflag && selector_.empty() && (msk.size() == 0 || msk.size() == ntrue) ) throw(AipsError("Trying to flag whole scantable.")); uChar userflag = 1 << 7; if ( unflag ) { userflag = 0 << 7; } if (whichrow > -1 ) { applyChanFlag(uInt(whichrow), msk, userflag); } else { for ( uInt i=0; i& msk, uChar flagval ) { if (whichrow >= table_.nrow() ) { throw( casa::indexError( whichrow, "asap::Scantable::applyChanFlag: Invalid row number" ) ); } Vector flgs = flagsCol_(whichrow); if ( msk.size() == 0 ) { flgs = flagval; flagsCol_.put(whichrow, flgs); return; } if ( int(msk.size()) != nchan() ) { throw(AipsError("Mask has incorrect number of channels.")); } if ( flgs.nelements() != msk.size() ) { throw(AipsError("Mask has incorrect number of channels." " Probably varying with IF. Please flag per IF")); } std::vector::const_iterator it; uInt j = 0; for (it = msk.begin(); it != msk.end(); ++it) { if ( *it ) { flgs(j) = flagval; } ++j; } flagsCol_.put(whichrow, flgs); } void Scantable::flagRow(const std::vector& rows, bool unflag) { if ( selector_.empty() && (rows.size() == table_.nrow()) ) throw(AipsError("Trying to flag whole scantable.")); uInt rowflag = (unflag ? 0 : 1); std::vector::const_iterator it; for (it = rows.begin(); it != rows.end(); ++it) flagrowCol_.put(*it, rowflag); } std::vector Scantable::getMask(int whichrow) const { Vector flags; flagsCol_.get(uInt(whichrow), flags); Vector bflag(flags.shape()); convertArray(bflag, flags); bflag = !bflag; std::vector mask; bflag.tovector(mask); return mask; } std::vector Scantable::getSpectrum( int whichrow, const std::string& poltype ) const { String ptype = poltype; if (poltype == "" ) ptype = getPolType(); if ( whichrow < 0 || whichrow >= nrow() ) throw(AipsError("Illegal row number.")); std::vector out; Vector arr; uInt requestedpol = polCol_(whichrow); String basetype = getPolType(); if ( ptype == basetype ) { specCol_.get(whichrow, arr); } else { CountedPtr stpol(STPol::getPolClass(Scantable::factories_, basetype)); uInt row = uInt(whichrow); stpol->setSpectra(getPolMatrix(row)); Float fang,fhand; fang = focusTable_.getTotalAngle(mfocusidCol_(row)); fhand = focusTable_.getFeedHand(mfocusidCol_(row)); stpol->setPhaseCorrections(fang, fhand); arr = stpol->getSpectrum(requestedpol, ptype); } if ( arr.nelements() == 0 ) pushLog("Not enough polarisations present to do the conversion."); arr.tovector(out); return out; } void Scantable::setSpectrum( const std::vector& spec, int whichrow ) { Vector spectrum(spec); Vector arr; specCol_.get(whichrow, arr); if ( spectrum.nelements() != arr.nelements() ) throw AipsError("The spectrum has incorrect number of channels."); specCol_.put(whichrow, spectrum); } String Scantable::generateName() { return (File::newUniqueName("./","temp")).baseName(); } const casa::Table& Scantable::table( ) const { return table_; } casa::Table& Scantable::table( ) { return table_; } std::string Scantable::getPolType() const { return table_.keywordSet().asString("POLTYPE"); } void Scantable::unsetSelection() { table_ = originalTable_; attach(); selector_.reset(); } void Scantable::setSelection( const STSelector& selection ) { Table tab = const_cast(selection).apply(originalTable_); if ( tab.nrow() == 0 ) { throw(AipsError("Selection contains no data. Not applying it.")); } table_ = tab; attach(); selector_ = selection; } std::string Scantable::summary( bool verbose ) { // Format header info ostringstream oss; oss << endl; oss << asap::SEPERATOR << endl; oss << " Scan Table Summary" << endl; oss << asap::SEPERATOR << endl; oss.flags(std::ios_base::left); oss << setw(15) << "Beams:" << setw(4) << nbeam() << endl << setw(15) << "IFs:" << setw(4) << nif() << endl << setw(15) << "Polarisations:" << setw(4) << npol() << "(" << getPolType() << ")" << endl << setw(15) << "Channels:" << nchan() << endl; String tmp; oss << setw(15) << "Observer:" << table_.keywordSet().asString("Observer") << endl; oss << setw(15) << "Obs Date:" << getTime(-1,true) << endl; table_.keywordSet().get("Project", tmp); oss << setw(15) << "Project:" << tmp << endl; table_.keywordSet().get("Obstype", tmp); oss << setw(15) << "Obs. Type:" << tmp << endl; table_.keywordSet().get("AntennaName", tmp); oss << setw(15) << "Antenna Name:" << tmp << endl; table_.keywordSet().get("FluxUnit", tmp); oss << setw(15) << "Flux Unit:" << tmp << endl; //Vector vec(moleculeTable_.getRestFrequencies()); int nid = moleculeTable_.nrow(); Bool firstline = True; oss << setw(15) << "Rest Freqs:"; for (int i=0; i 0) { Vector vec(moleculeTable_.getRestFrequency(i)); if (vec.nelements() > 0) { if (firstline) { oss << setprecision(10) << vec << " [Hz]" << endl; firstline=False; } else{ oss << setw(15)<<" " << setprecision(10) << vec << " [Hz]" << endl; } } else { oss << "none" << endl; } } } oss << setw(15) << "Abcissa:" << getAbcissaLabel(0) << endl; oss << selector_.print() << endl; oss << endl; // main table String dirtype = "Position (" + getDirectionRefString() + ")"; oss << setw(5) << "Scan" << setw(15) << "Source" << setw(10) << "Time" << setw(18) << "Integration" << setw(15) << "Source Type" << endl; oss << setw(5) << "" << setw(5) << "Beam" << setw(3) << "" << dirtype << endl; oss << setw(10) << "" << setw(3) << "IF" << setw(3) << "" << setw(8) << "Frame" << setw(16) << "RefVal" << setw(10) << "RefPix" << setw(12) << "Increment" << setw(7) << "Channels" << endl; oss << asap::SEPERATOR << endl; TableIterator iter(table_, "SCANNO"); while (!iter.pastEnd()) { Table subt = iter.table(); ROTableRow row(subt); MEpoch::ROScalarColumn timeCol(subt,"TIME"); const TableRecord& rec = row.get(0); oss << setw(4) << std::right << rec.asuInt("SCANNO") << std::left << setw(1) << "" << setw(15) << rec.asString("SRCNAME") << setw(10) << formatTime(timeCol(0), false); // count the cycles in the scan TableIterator cyciter(subt, "CYCLENO"); int nint = 0; while (!cyciter.pastEnd()) { ++nint; ++cyciter; } oss << setw(3) << std::right << nint << setw(3) << " x " << std::left << setw(11) << formatSec(rec.asFloat("INTERVAL")) << setw(1) << "" << setw(15) << SrcType::getName(rec.asInt("SRCTYPE")) << endl; TableIterator biter(subt, "BEAMNO"); while (!biter.pastEnd()) { Table bsubt = biter.table(); ROTableRow brow(bsubt); const TableRecord& brec = brow.get(0); uInt row0 = bsubt.rowNumbers(table_)[0]; oss << setw(5) << "" << setw(4) << std::right << brec.asuInt("BEAMNO")<< std::left; oss << setw(4) << "" << formatDirection(getDirection(row0)) << endl; TableIterator iiter(bsubt, "IFNO"); while (!iiter.pastEnd()) { Table isubt = iiter.table(); ROTableRow irow(isubt); const TableRecord& irec = irow.get(0); oss << setw(9) << ""; oss << setw(3) << std::right << irec.asuInt("IFNO") << std::left << setw(1) << "" << frequencies().print(irec.asuInt("FREQ_ID")) << setw(3) << "" << nchan(irec.asuInt("IFNO")) << endl; ++iiter; } ++biter; } ++iter; } /// @todo implement verbose mode return String(oss); } // std::string Scantable::getTime(int whichrow, bool showdate) const // { // MEpoch::ROScalarColumn timeCol(table_, "TIME"); // MEpoch me; // if (whichrow > -1) { // me = timeCol(uInt(whichrow)); // } else { // Double tm; // table_.keywordSet().get("UTC",tm); // me = MEpoch(MVEpoch(tm)); // } // return formatTime(me, showdate); // } std::string Scantable::getTime(int whichrow, bool showdate, uInt prec) const { MEpoch me; me = getEpoch(whichrow); return formatTime(me, showdate, prec); } MEpoch Scantable::getEpoch(int whichrow) const { if (whichrow > -1) { return timeCol_(uInt(whichrow)); } else { Double tm; table_.keywordSet().get("UTC",tm); return MEpoch(MVEpoch(tm)); } } std::string Scantable::getDirectionString(int whichrow) const { return formatDirection(getDirection(uInt(whichrow))); } SpectralCoordinate Scantable::getSpectralCoordinate(int whichrow) const { const MPosition& mp = getAntennaPosition(); const MDirection& md = getDirection(whichrow); const MEpoch& me = timeCol_(whichrow); //Double rf = moleculeTable_.getRestFrequency(mmolidCol_(whichrow)); Vector rf = moleculeTable_.getRestFrequency(mmolidCol_(whichrow)); return freqTable_.getSpectralCoordinate(md, mp, me, rf, mfreqidCol_(whichrow)); } std::vector< double > Scantable::getAbcissa( int whichrow ) const { if ( whichrow > int(table_.nrow()) ) throw(AipsError("Illegal row number")); std::vector stlout; int nchan = specCol_(whichrow).nelements(); String us = freqTable_.getUnitString(); if ( us == "" || us == "pixel" || us == "channel" ) { for (int i=0; i pixel(nchan); Vector world; indgen(pixel); if ( Unit(us) == Unit("Hz") ) { for ( int i=0; i < nchan; ++i) { Double world; spc.toWorld(world, pixel[i]); stlout.push_back(double(world)); } } else if ( Unit(us) == Unit("km/s") ) { Vector world; spc.pixelToVelocity(world, pixel); world.tovector(stlout); } return stlout; } void Scantable::setDirectionRefString( const std::string & refstr ) { MDirection::Types mdt; if (refstr != "" && !MDirection::getType(mdt, refstr)) { throw(AipsError("Illegal Direction frame.")); } if ( refstr == "" ) { String defaultstr = MDirection::showType(dirCol_.getMeasRef().getType()); table_.rwKeywordSet().define("DIRECTIONREF", defaultstr); } else { table_.rwKeywordSet().define("DIRECTIONREF", String(refstr)); } } std::string Scantable::getDirectionRefString( ) const { return table_.keywordSet().asString("DIRECTIONREF"); } MDirection Scantable::getDirection(int whichrow ) const { String usertype = table_.keywordSet().asString("DIRECTIONREF"); String type = MDirection::showType(dirCol_.getMeasRef().getType()); if ( usertype != type ) { MDirection::Types mdt; if (!MDirection::getType(mdt, usertype)) { throw(AipsError("Illegal Direction frame.")); } return dirCol_.convert(uInt(whichrow), mdt); } else { return dirCol_(uInt(whichrow)); } } std::string Scantable::getAbcissaLabel( int whichrow ) const { if ( whichrow > int(table_.nrow()) ) throw(AipsError("Illegal ro number")); const MPosition& mp = getAntennaPosition(); const MDirection& md = getDirection(whichrow); const MEpoch& me = timeCol_(whichrow); //const Double& rf = mmolidCol_(whichrow); const Vector rf = moleculeTable_.getRestFrequency(mmolidCol_(whichrow)); SpectralCoordinate spc = freqTable_.getSpectralCoordinate(md, mp, me, rf, mfreqidCol_(whichrow)); String s = "Channel"; Unit u = Unit(freqTable_.getUnitString()); if (u == Unit("km/s")) { s = CoordinateUtil::axisLabel(spc, 0, True,True, True); } else if (u == Unit("Hz")) { Vector wau(1);wau = u.getName(); spc.setWorldAxisUnits(wau); s = CoordinateUtil::axisLabel(spc, 0, True, True, False); } return s; } /** void asap::Scantable::setRestFrequencies( double rf, const std::string& name, const std::string& unit ) **/ void Scantable::setRestFrequencies( vector rf, const vector& name, const std::string& unit ) { ///@todo lookup in line table to fill in name and formattedname Unit u(unit); //Quantum urf(rf, u); Quantum >urf(rf, u); Vector formattedname(0); //cerr<<"Scantable::setRestFrequnecies="< tabvec(table_, "MOLECULE_ID"); tabvec = id; } /** void asap::Scantable::setRestFrequencies( const std::string& name ) { throw(AipsError("setRestFrequencies( const std::string& name ) NYI")); ///@todo implement } **/ void Scantable::setRestFrequencies( const vector& name ) { throw(AipsError("setRestFrequencies( const vector& name ) NYI")); ///@todo implement } std::vector< unsigned int > Scantable::rownumbers( ) const { std::vector stlout; Vector vec = table_.rowNumbers(); vec.tovector(stlout); return stlout; } Matrix Scantable::getPolMatrix( uInt whichrow ) const { ROTableRow row(table_); const TableRecord& rec = row.get(whichrow); Table t = originalTable_( originalTable_.col("SCANNO") == Int(rec.asuInt("SCANNO")) && originalTable_.col("BEAMNO") == Int(rec.asuInt("BEAMNO")) && originalTable_.col("IFNO") == Int(rec.asuInt("IFNO")) && originalTable_.col("CYCLENO") == Int(rec.asuInt("CYCLENO")) ); ROArrayColumn speccol(t, "SPECTRA"); return speccol.getColumn(); } std::vector< std::string > Scantable::columnNames( ) const { Vector vec = table_.tableDesc().columnNames(); return mathutil::tovectorstring(vec); } MEpoch::Types Scantable::getTimeReference( ) const { return MEpoch::castType(timeCol_.getMeasRef().getType()); } void Scantable::addFit( const STFitEntry& fit, int row ) { //cout << mfitidCol_(uInt(row)) << endl; LogIO os( LogOrigin( "Scantable", "addFit()", WHERE ) ) ; os << mfitidCol_(uInt(row)) << LogIO::POST ; uInt id = fitTable_.addEntry(fit, mfitidCol_(uInt(row))); mfitidCol_.put(uInt(row), id); } void Scantable::shift(int npix) { Vector fids(mfreqidCol_.getColumn()); genSort( fids, Sort::Ascending, Sort::QuickSort|Sort::NoDuplicates ); for (uInt i=0; i& scanlist) const { String tbpath; int ret = 0; if ( table_.keywordSet().isDefined("GBT_GO") ) { table_.keywordSet().get("GBT_GO", tbpath); Table t(tbpath,Table::Old); // check each scan if other scan of the pair exist int nscan = scanlist.size(); for (int i = 0; i < nscan; i++) { Table subt = t( t.col("SCAN") == scanlist[i]+1 ); if (subt.nrow()==0) { //cerr <<"Scan "<= nscan ) break; } } else { LogIO os( LogOrigin( "Scantable", "checkScanInfo()", WHERE ) ) ; //cerr<<"No reference to GBT_GO table."< Scantable::getDirectionVector(int whichrow) const { Vector Dir = dirCol_(whichrow).getAngle("rad").getValue(); std::vector dir; Dir.tovector(dir); return dir; } void asap::Scantable::reshapeSpectrum( int nmin, int nmax ) throw( casa::AipsError ) { // assumed that all rows have same nChan Vector arr = specCol_( 0 ) ; int nChan = arr.nelements() ; // if nmin < 0 or nmax < 0, nothing to do if ( nmin < 0 ) { throw( casa::indexError( nmin, "asap::Scantable::reshapeSpectrum: Invalid range. Negative index is specified." ) ) ; } if ( nmax < 0 ) { throw( casa::indexError( nmax, "asap::Scantable::reshapeSpectrum: Invalid range. Negative index is specified." ) ) ; } // if nmin > nmax, exchange values if ( nmin > nmax ) { int tmp = nmax ; nmax = nmin ; nmin = tmp ; LogIO os( LogOrigin( "Scantable", "reshapeSpectrum()", WHERE ) ) ; os << "Swap values. Applied range is [" << nmin << ", " << nmax << "]" << LogIO::POST ; } // if nmin exceeds nChan, nothing to do if ( nmin >= nChan ) { throw( casa::indexError( nmin, "asap::Scantable::reshapeSpectrum: Invalid range. Specified minimum exceeds nChan." ) ) ; } // if nmax exceeds nChan, reset nmax to nChan if ( nmax >= nChan ) { if ( nmin == 0 ) { // nothing to do LogIO os( LogOrigin( "Scantable", "reshapeSpectrum()", WHERE ) ) ; os << "Whole range is selected. Nothing to do." << LogIO::POST ; return ; } else { LogIO os( LogOrigin( "Scantable", "reshapeSpectrum()", WHERE ) ) ; os << "Specified maximum exceeds nChan. Applied range is [" << nmin << ", " << nChan-1 << "]." << LogIO::POST ; nmax = nChan - 1 ; } } // reshape specCol_ and flagCol_ for ( int irow = 0 ; irow < nrow() ; irow++ ) { reshapeSpectrum( nmin, nmax, irow ) ; } // update FREQUENCIES subtable Double refpix ; Double refval ; Double increment ; int freqnrow = freqTable_.table().nrow() ; Vector oldId( freqnrow ) ; Vector newId( freqnrow ) ; for ( int irow = 0 ; irow < freqnrow ; irow++ ) { freqTable_.getEntry( refpix, refval, increment, irow ) ; /*** * need to shift refpix to nmin * note that channel nmin in old index will be channel 0 in new one ***/ refval = refval - ( refpix - nmin ) * increment ; refpix = 0 ; freqTable_.setEntry( refpix, refval, increment, irow ) ; } // update nchan int newsize = nmax - nmin + 1 ; table_.rwKeywordSet().define( "nChan", newsize ) ; // update bandwidth // assumed all spectra in the scantable have same bandwidth table_.rwKeywordSet().define( "Bandwidth", increment * newsize ) ; return ; } void asap::Scantable::reshapeSpectrum( int nmin, int nmax, int irow ) { // reshape specCol_ and flagCol_ Vector oldspec = specCol_( irow ) ; Vector oldflag = flagsCol_( irow ) ; uInt newsize = nmax - nmin + 1 ; specCol_.put( irow, oldspec( Slice( nmin, newsize, 1 ) ) ) ; flagsCol_.put( irow, oldflag( Slice( nmin, newsize, 1 ) ) ) ; return ; } void asap::Scantable::regridChannel( int nChan, double dnu ) { LogIO os( LogOrigin( "Scantable", "regridChannel()", WHERE ) ) ; os << "Regrid abcissa with channel number " << nChan << " and spectral resoultion " << dnu << "Hz." << LogIO::POST ; // assumed that all rows have same nChan Vector arr = specCol_( 0 ) ; int oldsize = arr.nelements() ; // if oldsize == nChan, nothing to do if ( oldsize == nChan ) { os << "Specified channel number is same as current one. Nothing to do." << LogIO::POST ; return ; } // if oldChan < nChan, unphysical operation if ( oldsize < nChan ) { os << "Unphysical operation. Nothing to do." << LogIO::POST ; return ; } // change channel number for specCol_ and flagCol_ Vector newspec( nChan, 0 ) ; Vector newflag( nChan, false ) ; vector coordinfo = getCoordInfo() ; string oldinfo = coordinfo[0] ; coordinfo[0] = "Hz" ; setCoordInfo( coordinfo ) ; for ( int irow = 0 ; irow < nrow() ; irow++ ) { regridChannel( nChan, dnu, irow ) ; } coordinfo[0] = oldinfo ; setCoordInfo( coordinfo ) ; // NOTE: this method does not update metadata such as // FREQUENCIES subtable, nChan, Bandwidth, etc. return ; } void asap::Scantable::regridChannel( int nChan, double dnu, int irow ) { // logging //ofstream ofs( "average.log", std::ios::out | std::ios::app ) ; //ofs << "IFNO = " << getIF( irow ) << " irow = " << irow << endl ; Vector oldspec = specCol_( irow ) ; Vector oldflag = flagsCol_( irow ) ; Vector newspec( nChan, 0 ) ; Vector newflag( nChan, false ) ; // regrid vector abcissa = getAbcissa( irow ) ; int oldsize = abcissa.size() ; double olddnu = abcissa[1] - abcissa[0] ; //int refChan = 0 ; //double frac = 0.0 ; //double wedge = 0.0 ; //double pile = 0.0 ; int ichan = 0 ; double wsum = 0.0 ; Vector z( nChan ) ; z[0] = abcissa[0] - 0.5 * olddnu + 0.5 * dnu ; for ( int ii = 1 ; ii < nChan ; ii++ ) z[ii] = z[ii-1] + dnu ; Vector zi( nChan+1 ) ; Vector yi( oldsize + 1 ) ; zi[0] = z[0] - 0.5 * dnu ; zi[1] = z[0] + 0.5 * dnu ; for ( int ii = 2 ; ii < nChan ; ii++ ) zi[ii] = zi[ii-1] + dnu ; zi[nChan] = z[nChan-1] + 0.5 * dnu ; yi[0] = abcissa[0] - 0.5 * olddnu ; yi[1] = abcissa[1] + 0.5 * olddnu ; for ( int ii = 2 ; ii < oldsize ; ii++ ) yi[ii] = abcissa[ii-1] + olddnu ; yi[oldsize] = abcissa[oldsize-1] + 0.5 * olddnu ; if ( dnu > 0.0 ) { for ( int ii = 0 ; ii < nChan ; ii++ ) { double zl = zi[ii] ; double zr = zi[ii+1] ; for ( int j = ichan ; j < oldsize ; j++ ) { double yl = yi[j] ; double yr = yi[j+1] ; if ( yl <= zl ) { if ( yr <= zl ) { continue ; } else if ( yr <= zr ) { newspec[ii] += oldspec[j] * ( yr - zl ) ; newflag[ii] = newflag[ii] || oldflag[j] ; wsum += ( yr - zl ) ; } else { newspec[ii] += oldspec[j] * dnu ; newflag[ii] = newflag[ii] || oldflag[j] ; wsum += dnu ; ichan = j ; break ; } } else if ( yl < zr ) { if ( yr <= zr ) { newspec[ii] += oldspec[j] * ( yr - yl ) ; newflag[ii] = newflag[ii] || oldflag[j] ; wsum += ( yr - yl ) ; } else { newspec[ii] += oldspec[j] * ( zr - yl ) ; newflag[ii] = newflag[ii] || oldflag[j] ; wsum += ( zr - yl ) ; ichan = j ; break ; } } else { ichan = j - 1 ; break ; } } newspec[ii] /= wsum ; wsum = 0.0 ; } } else if ( dnu < 0.0 ) { for ( int ii = 0 ; ii < nChan ; ii++ ) { double zl = zi[ii] ; double zr = zi[ii+1] ; for ( int j = ichan ; j < oldsize ; j++ ) { double yl = yi[j] ; double yr = yi[j+1] ; if ( yl >= zl ) { if ( yr >= zl ) { continue ; } else if ( yr >= zr ) { newspec[ii] += oldspec[j] * abs( yr - zl ) ; newflag[ii] = newflag[ii] || oldflag[j] ; wsum += abs( yr - zl ) ; } else { newspec[ii] += oldspec[j] * abs( dnu ) ; newflag[ii] = newflag[ii] || oldflag[j] ; wsum += abs( dnu ) ; ichan = j ; break ; } } else if ( yl > zr ) { if ( yr >= zr ) { newspec[ii] += oldspec[j] * abs( yr - yl ) ; newflag[ii] = newflag[ii] || oldflag[j] ; wsum += abs( yr - yl ) ; } else { newspec[ii] += oldspec[j] * abs( zr - yl ) ; newflag[ii] = newflag[ii] || oldflag[j] ; wsum += abs( zr - yl ) ; ichan = j ; break ; } } else { ichan = j - 1 ; break ; } } newspec[ii] /= wsum ; wsum = 0.0 ; } } // * ichan = 0 // ***/ // //ofs << "olddnu = " << olddnu << ", dnu = " << dnu << endl ; // pile += dnu ; // wedge = olddnu * ( refChan + 1 ) ; // while ( wedge < pile ) { // newspec[0] += olddnu * oldspec[refChan] ; // newflag[0] = newflag[0] || oldflag[refChan] ; // //ofs << "channel " << refChan << " is included in new channel 0" << endl ; // refChan++ ; // wedge += olddnu ; // wsum += olddnu ; // //ofs << "newspec[0] = " << newspec[0] << " wsum = " << wsum << endl ; // } // frac = ( wedge - pile ) / olddnu ; // wsum += ( 1.0 - frac ) * olddnu ; // newspec[0] += ( 1.0 - frac ) * olddnu * oldspec[refChan] ; // newflag[0] = newflag[0] || oldflag[refChan] ; // //ofs << "channel " << refChan << " is partly included in new channel 0" << " with fraction of " << ( 1.0 - frac ) << endl ; // //ofs << "newspec[0] = " << newspec[0] << " wsum = " << wsum << endl ; // newspec[0] /= wsum ; // //ofs << "newspec[0] = " << newspec[0] << endl ; // //ofs << "wedge = " << wedge << ", pile = " << pile << endl ; // /*** // * ichan = 1 - nChan-2 // ***/ // for ( int ichan = 1 ; ichan < nChan - 1 ; ichan++ ) { // pile += dnu ; // newspec[ichan] += frac * olddnu * oldspec[refChan] ; // newflag[ichan] = newflag[ichan] || oldflag[refChan] ; // //ofs << "channel " << refChan << " is partly included in new channel " << ichan << " with fraction of " << frac << endl ; // refChan++ ; // wedge += olddnu ; // wsum = frac * olddnu ; // //ofs << "newspec[" << ichan << "] = " << newspec[ichan] << " wsum = " << wsum << endl ; // while ( wedge < pile ) { // newspec[ichan] += olddnu * oldspec[refChan] ; // newflag[ichan] = newflag[ichan] || oldflag[refChan] ; // //ofs << "channel " << refChan << " is included in new channel " << ichan << endl ; // refChan++ ; // wedge += olddnu ; // wsum += olddnu ; // //ofs << "newspec[" << ichan << "] = " << newspec[ichan] << " wsum = " << wsum << endl ; // } // frac = ( wedge - pile ) / olddnu ; // wsum += ( 1.0 - frac ) * olddnu ; // newspec[ichan] += ( 1.0 - frac ) * olddnu * oldspec[refChan] ; // newflag[ichan] = newflag[ichan] || oldflag[refChan] ; // //ofs << "channel " << refChan << " is partly included in new channel " << ichan << " with fraction of " << ( 1.0 - frac ) << endl ; // //ofs << "wedge = " << wedge << ", pile = " << pile << endl ; // //ofs << "newspec[" << ichan << "] = " << newspec[ichan] << " wsum = " << wsum << endl ; // newspec[ichan] /= wsum ; // //ofs << "newspec[" << ichan << "] = " << newspec[ichan] << endl ; // } // /*** // * ichan = nChan-1 // ***/ // // NOTE: Assumed that all spectra have the same bandwidth // pile += dnu ; // newspec[nChan-1] += frac * olddnu * oldspec[refChan] ; // newflag[nChan-1] = newflag[nChan-1] || oldflag[refChan] ; // //ofs << "channel " << refChan << " is partly included in new channel " << nChan-1 << " with fraction of " << frac << endl ; // refChan++ ; // wedge += olddnu ; // wsum = frac * olddnu ; // //ofs << "newspec[" << nChan - 1 << "] = " << newspec[nChan-1] << " wsum = " << wsum << endl ; // for ( int jchan = refChan ; jchan < oldsize ; jchan++ ) { // newspec[nChan-1] += olddnu * oldspec[jchan] ; // newflag[nChan-1] = newflag[nChan-1] || oldflag[jchan] ; // wsum += olddnu ; // //ofs << "channel " << jchan << " is included in new channel " << nChan-1 << " with fraction of " << frac << endl ; // //ofs << "newspec[" << nChan - 1 << "] = " << newspec[nChan-1] << " wsum = " << wsum << endl ; // } // //ofs << "wedge = " << wedge << ", pile = " << pile << endl ; // //ofs << "newspec[" << nChan - 1 << "] = " << newspec[nChan-1] << " wsum = " << wsum << endl ; // newspec[nChan-1] /= wsum ; // //ofs << "newspec[" << nChan - 1 << "] = " << newspec[nChan-1] << endl ; // // ofs.close() ; specCol_.put( irow, newspec ) ; flagsCol_.put( irow, newflag ) ; return ; } std::vector Scantable::getWeather(int whichrow) const { std::vector out(5); //Float temperature, pressure, humidity, windspeed, windaz; weatherTable_.getEntry(out[0], out[1], out[2], out[3], out[4], mweatheridCol_(uInt(whichrow))); return out; } bool Scantable::getFlagtraFast(uInt whichrow) { uChar flag; Vector flags; flagsCol_.get(whichrow, flags); flag = flags[0]; for (uInt i = 1; i < flags.size(); ++i) { flag &= flags[i]; } return ((flag >> 7) == 1); } void Scantable::polyBaseline(const std::vector& mask, int order, bool outLogger, const std::string& blfile) { ofstream ofs; String coordInfo; bool hasSameNchan = true; bool outTextFile = false; if (blfile != "") { ofs.open(blfile.c_str(), ios::out | ios::app); if (ofs) outTextFile = true; } if (outLogger || outTextFile) { coordInfo = getCoordInfo()[0]; if (coordInfo == "") coordInfo = "channel"; hasSameNchan = hasSameNchanOverIFs(); } Fitter fitter = Fitter(); fitter.setExpression("poly", order); int nRow = nrow(); std::vector chanMask; for (int whichrow = 0; whichrow < nRow; ++whichrow) { chanMask = getCompositeChanMask(whichrow, mask); fitBaseline(chanMask, whichrow, fitter); setSpectrum(fitter.getResidual(), whichrow); outputFittingResult(outLogger, outTextFile, chanMask, whichrow, coordInfo, hasSameNchan, ofs, "polyBaseline()", fitter); } if (outTextFile) ofs.close(); } void Scantable::autoPolyBaseline(const std::vector& mask, int order, const std::vector& edge, float threshold, int chanAvgLimit, bool outLogger, const std::string& blfile) { ofstream ofs; String coordInfo; bool hasSameNchan = true; bool outTextFile = false; if (blfile != "") { ofs.open(blfile.c_str(), ios::out | ios::app); if (ofs) outTextFile = true; } if (outLogger || outTextFile) { coordInfo = getCoordInfo()[0]; if (coordInfo == "") coordInfo = "channel"; hasSameNchan = hasSameNchanOverIFs(); } Fitter fitter = Fitter(); fitter.setExpression("poly", order); int nRow = nrow(); std::vector chanMask; int minEdgeSize = getIFNos().size()*2; STLineFinder lineFinder = STLineFinder(); lineFinder.setOptions(threshold, 3, chanAvgLimit); for (int whichrow = 0; whichrow < nRow; ++whichrow) { //------------------------------------------------------- //chanMask = getCompositeChanMask(whichrow, mask, edge, minEdgeSize, lineFinder); //------------------------------------------------------- int edgeSize = edge.size(); std::vector currentEdge; if (edgeSize >= 2) { int idx = 0; if (edgeSize > 2) { if (edgeSize < minEdgeSize) { throw(AipsError("Length of edge element info is less than that of IFs")); } idx = 2 * getIF(whichrow); } currentEdge.push_back(edge[idx]); currentEdge.push_back(edge[idx+1]); } else { throw(AipsError("Wrong length of edge element")); } lineFinder.setData(getSpectrum(whichrow)); lineFinder.findLines(getCompositeChanMask(whichrow, mask), currentEdge, whichrow); chanMask = lineFinder.getMask(); //------------------------------------------------------- fitBaseline(chanMask, whichrow, fitter); setSpectrum(fitter.getResidual(), whichrow); outputFittingResult(outLogger, outTextFile, chanMask, whichrow, coordInfo, hasSameNchan, ofs, "autoPolyBaseline()", fitter); } if (outTextFile) ofs.close(); } void Scantable::cubicSplineBaseline(const std::vector& mask, int nPiece, float thresClip, int nIterClip, bool outLogger, const std::string& blfile) { ofstream ofs; String coordInfo; bool hasSameNchan = true; bool outTextFile = false; if (blfile != "") { ofs.open(blfile.c_str(), ios::out | ios::app); if (ofs) outTextFile = true; } if (outLogger || outTextFile) { coordInfo = getCoordInfo()[0]; if (coordInfo == "") coordInfo = "channel"; hasSameNchan = hasSameNchanOverIFs(); } //Fitter fitter = Fitter(); //fitter.setExpression("cspline", nPiece); int nRow = nrow(); std::vector chanMask; for (int whichrow = 0; whichrow < nRow; ++whichrow) { chanMask = getCompositeChanMask(whichrow, mask); //fitBaseline(chanMask, whichrow, fitter, thresClip, nIterClip); //setSpectrum(fitter.getResidual(), whichrow); std::vector pieceRanges; std::vector params; std::vector res = doCubicSplineFitting(getSpectrum(whichrow), chanMask, pieceRanges, params, nPiece, thresClip, nIterClip, true); setSpectrum(res, whichrow); // std::vector fixed; fixed.clear(); outputFittingResult(outLogger, outTextFile, chanMask, whichrow, coordInfo, hasSameNchan, ofs, "cubicSplineBaseline()", pieceRanges, params, fixed); } if (outTextFile) ofs.close(); } void Scantable::autoCubicSplineBaseline(const std::vector& mask, int nPiece, float thresClip, int nIterClip, const std::vector& edge, float threshold, int chanAvgLimit, bool outLogger, const std::string& blfile) { ofstream ofs; String coordInfo; bool hasSameNchan = true; bool outTextFile = false; if (blfile != "") { ofs.open(blfile.c_str(), ios::out | ios::app); if (ofs) outTextFile = true; } if (outLogger || outTextFile) { coordInfo = getCoordInfo()[0]; if (coordInfo == "") coordInfo = "channel"; hasSameNchan = hasSameNchanOverIFs(); } //Fitter fitter = Fitter(); //fitter.setExpression("cspline", nPiece); int nRow = nrow(); std::vector chanMask; int minEdgeSize = getIFNos().size()*2; STLineFinder lineFinder = STLineFinder(); lineFinder.setOptions(threshold, 3, chanAvgLimit); for (int whichrow = 0; whichrow < nRow; ++whichrow) { //------------------------------------------------------- //chanMask = getCompositeChanMask(whichrow, mask, edge, minEdgeSize, lineFinder); //------------------------------------------------------- int edgeSize = edge.size(); std::vector currentEdge; if (edgeSize >= 2) { int idx = 0; if (edgeSize > 2) { if (edgeSize < minEdgeSize) { throw(AipsError("Length of edge element info is less than that of IFs")); } idx = 2 * getIF(whichrow); } currentEdge.push_back(edge[idx]); currentEdge.push_back(edge[idx+1]); } else { throw(AipsError("Wrong length of edge element")); } lineFinder.setData(getSpectrum(whichrow)); lineFinder.findLines(getCompositeChanMask(whichrow, mask), currentEdge, whichrow); chanMask = lineFinder.getMask(); //------------------------------------------------------- //fitBaseline(chanMask, whichrow, fitter, thresClip, nIterClip); //setSpectrum(fitter.getResidual(), whichrow); std::vector pieceRanges; std::vector params; std::vector res = doCubicSplineFitting(getSpectrum(whichrow), chanMask, pieceRanges, params, nPiece, thresClip, nIterClip, true); setSpectrum(res, whichrow); // std::vector fixed; fixed.clear(); outputFittingResult(outLogger, outTextFile, chanMask, whichrow, coordInfo, hasSameNchan, ofs, "autoCubicSplineBaseline()", pieceRanges, params, fixed); } if (outTextFile) ofs.close(); } std::vector Scantable::doCubicSplineFitting(const std::vector& data, const std::vector& mask, std::vector& idxEdge, std::vector& params, int nPiece, float thresClip, int nIterClip, bool getResidual) { if (nPiece < 1) { throw(AipsError("wrong number of the sections for fitting")); } if (data.size() != mask.size()) { throw(AipsError("data and mask have different size")); } int nChan = data.size(); std::vector maskArray; std::vector x; for (int i = 0; i < nChan; ++i) { maskArray.push_back(mask[i] ? 1 : 0); if (mask[i]) { x.push_back(i); } } int nData = x.size(); int nElement = (int)(floor(floor((double)(nData/nPiece))+0.5)); std::vector invEdge; idxEdge.clear(); idxEdge.push_back(x[0]); for (int i = 1; i < nPiece; ++i) { int valX = x[nElement*i]; idxEdge.push_back(valX); invEdge.push_back(1.0/(double)valX); } idxEdge.push_back(x[x.size()-1]+1); std::vector x1, x2, x3, z1, x1z1, x2z1, x3z1, r1; for (int i = 0; i < nChan; ++i) { double di = (double)i; double dD = (double)data[i]; x1.push_back(di); x2.push_back(di*di); x3.push_back(di*di*di); z1.push_back(dD); x1z1.push_back(dD*di); x2z1.push_back(dD*di*di); x3z1.push_back(dD*di*di*di); r1.push_back(0.0); } int currentNData = nData; int nDOF = nPiece + 3; //number of parameters to solve, namely, 4+(nPiece-1). for (int nClip = 0; nClip < nIterClip+1; ++nClip) { // xMatrix : horizontal concatenation of // the least-sq. matrix (left) and an // identity matrix (right). // the right part is used to calculate the inverse matrix of the left part. double xMatrix[nDOF][2*nDOF]; double zMatrix[nDOF]; for (int i = 0; i < nDOF; ++i) { for (int j = 0; j < 2*nDOF; ++j) { xMatrix[i][j] = 0.0; } xMatrix[i][nDOF+i] = 1.0; zMatrix[i] = 0.0; } for (int n = 0; n < nPiece; ++n) { for (int i = idxEdge[n]; i < idxEdge[n+1]; ++i) { if (maskArray[i] == 0) continue; xMatrix[0][0] += 1.0; xMatrix[0][1] += x1[i]; xMatrix[0][2] += x2[i]; xMatrix[0][3] += x3[i]; xMatrix[1][1] += x2[i]; xMatrix[1][2] += x3[i]; xMatrix[1][3] += x2[i]*x2[i]; xMatrix[2][2] += x2[i]*x2[i]; xMatrix[2][3] += x3[i]*x2[i]; xMatrix[3][3] += x3[i]*x3[i]; zMatrix[0] += z1[i]; zMatrix[1] += x1z1[i]; zMatrix[2] += x2z1[i]; zMatrix[3] += x3z1[i]; for (int j = 0; j < n; ++j) { double q = 1.0 - x1[i]*invEdge[j]; q = q*q*q; xMatrix[0][j+4] += q; xMatrix[1][j+4] += q*x1[i]; xMatrix[2][j+4] += q*x2[i]; xMatrix[3][j+4] += q*x3[i]; for (int k = 0; k < j; ++k) { double r = 1.0 - x1[i]*invEdge[k]; r = r*r*r; xMatrix[k+4][j+4] += r*q; } xMatrix[j+4][j+4] += q*q; zMatrix[j+4] += q*z1[i]; } } } for (int i = 0; i < nDOF; ++i) { for (int j = 0; j < i; ++j) { xMatrix[i][j] = xMatrix[j][i]; } } std::vector invDiag; for (int i = 0; i < nDOF; ++i) { invDiag.push_back(1.0/xMatrix[i][i]); for (int j = 0; j < nDOF; ++j) { xMatrix[i][j] *= invDiag[i]; } } for (int k = 0; k < nDOF; ++k) { for (int i = 0; i < nDOF; ++i) { if (i != k) { double factor1 = xMatrix[k][k]; double factor2 = xMatrix[i][k]; for (int j = k; j < 2*nDOF; ++j) { xMatrix[i][j] *= factor1; xMatrix[i][j] -= xMatrix[k][j]*factor2; xMatrix[i][j] /= factor1; } } } double xDiag = xMatrix[k][k]; for (int j = k; j < 2*nDOF; ++j) { xMatrix[k][j] /= xDiag; } } for (int i = 0; i < nDOF; ++i) { for (int j = 0; j < nDOF; ++j) { xMatrix[i][nDOF+j] *= invDiag[j]; } } //compute a vector y which consists of the coefficients of the best-fit spline curves //(a0,a1,a2,a3(,b3,c3,...)), namely, the ones for the leftmost piece and the ones of //cubic terms for the other pieces (in case nPiece>1). std::vector y; y.clear(); for (int i = 0; i < nDOF; ++i) { y.push_back(0.0); for (int j = 0; j < nDOF; ++j) { y[i] += xMatrix[i][nDOF+j]*zMatrix[j]; } } double a0 = y[0]; double a1 = y[1]; double a2 = y[2]; double a3 = y[3]; params.clear(); for (int n = 0; n < nPiece; ++n) { for (int i = idxEdge[n]; i < idxEdge[n+1]; ++i) { r1[i] = a0 + a1*x1[i] + a2*x2[i] + a3*x3[i]; } params.push_back(a0); params.push_back(a1); params.push_back(a2); params.push_back(a3); if (n == nPiece-1) break; double d = y[4+n]; double iE = invEdge[n]; a0 += d; a1 -= 3.0*d*iE; a2 += 3.0*d*iE*iE; a3 -= d*iE*iE*iE; } if ((nClip == nIterClip) || (thresClip <= 0.0)) { break; } else { std::vector rz; double stdDev = 0.0; for (int i = 0; i < nChan; ++i) { double val = abs(r1[i] - z1[i]); rz.push_back(val); stdDev += val*val*(double)maskArray[i]; } stdDev = sqrt(stdDev/(double)nData); double thres = stdDev * thresClip; int newNData = 0; for (int i = 0; i < nChan; ++i) { if (rz[i] >= thres) { maskArray[i] = 0; } if (maskArray[i] > 0) { newNData++; } } if (newNData == currentNData) { break; //no more flag to add. iteration stops. } else { currentNData = newNData; } } } std::vector result; if (getResidual) { for (int i = 0; i < nChan; ++i) { result.push_back((float)(z1[i] - r1[i])); } } else { for (int i = 0; i < nChan; ++i) { result.push_back((float)r1[i]); } } return result; } void Scantable::sinusoidBaseline(const std::vector& mask, int nMinWavesInSW, int nMaxWavesInSW, float thresClip, int nIterClip, bool outLogger, const std::string& blfile) { ofstream ofs; String coordInfo; bool hasSameNchan = true; bool outTextFile = false; if (blfile != "") { ofs.open(blfile.c_str(), ios::out | ios::app); if (ofs) outTextFile = true; } if (outLogger || outTextFile) { coordInfo = getCoordInfo()[0]; if (coordInfo == "") coordInfo = "channel"; hasSameNchan = hasSameNchanOverIFs(); } //Fitter fitter = Fitter(); //fitter.setExpression("sinusoid", nMaxWavesInSW); int nRow = nrow(); std::vector chanMask; for (int whichrow = 0; whichrow < nRow; ++whichrow) { chanMask = getCompositeChanMask(whichrow, mask); //fitBaseline(chanMask, whichrow, fitter, thresClip, nIterClip); //setSpectrum(fitter.getResidual(), whichrow); std::vector pieceRanges; std::vector params; std::vector res = doSinusoidFitting(getSpectrum(whichrow), chanMask, nMinWavesInSW, nMaxWavesInSW, params, thresClip, nIterClip, true); setSpectrum(res, whichrow); // std::vector fixed; fixed.clear(); outputFittingResult(outLogger, outTextFile, chanMask, whichrow, coordInfo, hasSameNchan, ofs, "sinusoidBaseline()", pieceRanges, params, fixed); } if (outTextFile) ofs.close(); } void Scantable::autoSinusoidBaseline(const std::vector& mask, int nMinWavesInSW, int nMaxWavesInSW, float thresClip, int nIterClip, const std::vector& edge, float threshold, int chanAvgLimit, bool outLogger, const std::string& blfile) { ofstream ofs; String coordInfo; bool hasSameNchan = true; bool outTextFile = false; if (blfile != "") { ofs.open(blfile.c_str(), ios::out | ios::app); if (ofs) outTextFile = true; } if (outLogger || outTextFile) { coordInfo = getCoordInfo()[0]; if (coordInfo == "") coordInfo = "channel"; hasSameNchan = hasSameNchanOverIFs(); } //Fitter fitter = Fitter(); //fitter.setExpression("sinusoid", nMaxWavesInSW); int nRow = nrow(); std::vector chanMask; int minEdgeSize = getIFNos().size()*2; STLineFinder lineFinder = STLineFinder(); lineFinder.setOptions(threshold, 3, chanAvgLimit); for (int whichrow = 0; whichrow < nRow; ++whichrow) { //------------------------------------------------------- //chanMask = getCompositeChanMask(whichrow, mask, edge, minEdgeSize, lineFinder); //------------------------------------------------------- int edgeSize = edge.size(); std::vector currentEdge; if (edgeSize >= 2) { int idx = 0; if (edgeSize > 2) { if (edgeSize < minEdgeSize) { throw(AipsError("Length of edge element info is less than that of IFs")); } idx = 2 * getIF(whichrow); } currentEdge.push_back(edge[idx]); currentEdge.push_back(edge[idx+1]); } else { throw(AipsError("Wrong length of edge element")); } lineFinder.setData(getSpectrum(whichrow)); lineFinder.findLines(getCompositeChanMask(whichrow, mask), currentEdge, whichrow); chanMask = lineFinder.getMask(); //------------------------------------------------------- //fitBaseline(chanMask, whichrow, fitter, thresClip, nIterClip); //setSpectrum(fitter.getResidual(), whichrow); std::vector pieceRanges; std::vector params; std::vector res = doSinusoidFitting(getSpectrum(whichrow), chanMask, nMinWavesInSW, nMaxWavesInSW, params, thresClip, nIterClip, true); setSpectrum(res, whichrow); // std::vector fixed; fixed.clear(); outputFittingResult(outLogger, outTextFile, chanMask, whichrow, coordInfo, hasSameNchan, ofs, "autoSinusoidBaseline()", pieceRanges, params, fixed); } if (outTextFile) ofs.close(); } std::vector Scantable::doSinusoidFitting(const std::vector& data, const std::vector& mask, int nMinWavesInSW, int nMaxWavesInSW, std::vector& params, float thresClip, int nIterClip, bool getResidual) { if (data.size() != mask.size()) { throw(AipsError("data and mask have different size")); } if ((nMinWavesInSW < 0) || (nMaxWavesInSW < 0)) { throw(AipsError("wave number must be positive or zero (i.e. constant)")); } else { if (nMaxWavesInSW < nMinWavesInSW) { int temp = nMaxWavesInSW; nMaxWavesInSW = nMinWavesInSW; nMinWavesInSW = temp; } } int nChan = data.size(); std::vector maskArray; std::vector x; for (int i = 0; i < nChan; ++i) { maskArray.push_back(mask[i] ? 1 : 0); if (mask[i]) { x.push_back(i); } } int nData = x.size(); std::vector x1, x2, x3, z1, x1z1, x2z1, x3z1, r1; for (int i = 0; i < nChan; ++i) { double di = (double)i; double dD = (double)data[i]; x1.push_back(di); x2.push_back(di*di); x3.push_back(di*di*di); z1.push_back(dD); x1z1.push_back(di*dD); x2z1.push_back(di*di*dD); x3z1.push_back(di*di*di*dD); r1.push_back(0.0); } int currentNData = nData; int nDOF = nMaxWavesInSW + 1; //number of parameters to solve. for (int nClip = 0; nClip < nIterClip+1; ++nClip) { //xMatrix : horizontal concatenation of the least-sq. matrix (left) and an //identity matrix (right). //the right part is used to calculate the inverse matrix of the left part. double xMatrix[nDOF][2*nDOF]; double zMatrix[nDOF]; for (int i = 0; i < nDOF; ++i) { for (int j = 0; j < 2*nDOF; ++j) { xMatrix[i][j] = 0.0; } xMatrix[i][nDOF+i] = 1.0; zMatrix[i] = 0.0; } for (int i = 0; i < currentNData; ++i) { if (maskArray[i] == 0) continue; xMatrix[0][0] += 1.0; xMatrix[0][1] += x1[i]; xMatrix[0][2] += x2[i]; xMatrix[0][3] += x3[i]; xMatrix[1][1] += x2[i]; xMatrix[1][2] += x3[i]; xMatrix[1][3] += x2[i]*x2[i]; xMatrix[2][2] += x2[i]*x2[i]; xMatrix[2][3] += x3[i]*x2[i]; xMatrix[3][3] += x3[i]*x3[i]; zMatrix[0] += z1[i]; zMatrix[1] += x1z1[i]; zMatrix[2] += x2z1[i]; zMatrix[3] += x3z1[i]; } for (int i = 0; i < nDOF; ++i) { for (int j = 0; j < i; ++j) { xMatrix[i][j] = xMatrix[j][i]; } } std::vector invDiag; for (int i = 0; i < nDOF; ++i) { invDiag.push_back(1.0/xMatrix[i][i]); for (int j = 0; j < nDOF; ++j) { xMatrix[i][j] *= invDiag[i]; } } for (int k = 0; k < nDOF; ++k) { for (int i = 0; i < nDOF; ++i) { if (i != k) { double factor1 = xMatrix[k][k]; double factor2 = xMatrix[i][k]; for (int j = k; j < 2*nDOF; ++j) { xMatrix[i][j] *= factor1; xMatrix[i][j] -= xMatrix[k][j]*factor2; xMatrix[i][j] /= factor1; } } } double xDiag = xMatrix[k][k]; for (int j = k; j < 2*nDOF; ++j) { xMatrix[k][j] /= xDiag; } } for (int i = 0; i < nDOF; ++i) { for (int j = 0; j < nDOF; ++j) { xMatrix[i][nDOF+j] *= invDiag[j]; } } //compute a vector y which consists of the coefficients of the sinusoids forming the //best-fit curves (a0,b0,a1,b1,...), namely, a* and b* are of sine and cosine functions, //respectively. std::vector y; for (int i = 0; i < nDOF; ++i) { y.push_back(0.0); for (int j = 0; j < nDOF; ++j) { y[i] += xMatrix[i][nDOF+j]*zMatrix[j]; } } double a0 = y[0]; double a1 = y[1]; double a2 = y[2]; double a3 = y[3]; params.clear(); for (int i = 0; i < nChan; ++i) { r1[i] = a0 + a1*x1[i] + a2*x2[i] + a3*x3[i]; } params.push_back(a0); params.push_back(a1); params.push_back(a2); params.push_back(a3); if ((nClip == nIterClip) || (thresClip <= 0.0)) { break; } else { std::vector rz; double stdDev = 0.0; for (int i = 0; i < nChan; ++i) { double val = abs(r1[i] - z1[i]); rz.push_back(val); stdDev += val*val*(double)maskArray[i]; } stdDev = sqrt(stdDev/(double)nData); double thres = stdDev * thresClip; int newNData = 0; for (int i = 0; i < nChan; ++i) { if (rz[i] >= thres) { maskArray[i] = 0; } if (maskArray[i] > 0) { newNData++; } } if (newNData == currentNData) { break; //no additional flag. finish iteration. } else { currentNData = newNData; } } } std::vector result; if (getResidual) { for (int i = 0; i < nChan; ++i) { result.push_back((float)(z1[i] - r1[i])); } } else { for (int i = 0; i < nChan; ++i) { result.push_back((float)r1[i]); } } return result; } void Scantable::fitBaseline(const std::vector& mask, int whichrow, Fitter& fitter) { std::vector abcsd = getAbcissa(whichrow); std::vector abcs; for (uInt i = 0; i < abcsd.size(); ++i) { abcs.push_back((float)abcsd[i]); } std::vector spec = getSpectrum(whichrow); fitter.setData(abcs, spec, mask); fitter.lfit(); } std::vector Scantable::getCompositeChanMask(int whichrow, const std::vector& inMask) { std::vector chanMask = getMask(whichrow); uInt chanMaskSize = chanMask.size(); if (chanMaskSize != inMask.size()) { throw(AipsError("different mask sizes")); } for (uInt i = 0; i < chanMaskSize; ++i) { chanMask[i] = chanMask[i] && inMask[i]; } return chanMask; } /* std::vector Scantable::getCompositeChanMask(int whichrow, const std::vector& inMask, const std::vector& edge, const int minEdgeSize, STLineFinder& lineFinder) { int edgeSize = edge.size(); std::vector currentEdge; if (edgeSize >= 2) { int idx = 0; if (edgeSize > 2) { if (edgeSize < minEdgeSize) { throw(AipsError("Length of edge element info is less than that of IFs")); } idx = 2 * getIF(whichrow); } currentEdge.push_back(edge[idx]); currentEdge.push_back(edge[idx+1]); } else { throw(AipsError("Wrong length of edge element")); } lineFinder.setData(getSpectrum(whichrow)); lineFinder.findLines(getCompositeChanMask(whichrow, inMask), currentEdge, whichrow); return lineFinder.getMask(); } */ /* for poly. the variations of outputFittingResult() should be merged into one eventually (2011/3/10 WK) */ void Scantable::outputFittingResult(bool outLogger, bool outTextFile, const std::vector& chanMask, int whichrow, const casa::String& coordInfo, bool hasSameNchan, ofstream& ofs, const casa::String& funcName, Fitter& fitter) { if (outLogger || outTextFile) { std::vector params = fitter.getParameters(); std::vector fixed = fitter.getFixedParameters(); float rms = getRms(chanMask, whichrow); String masklist = getMaskRangeList(chanMask, whichrow, coordInfo, hasSameNchan); if (outLogger) { LogIO ols(LogOrigin("Scantable", funcName, WHERE)); ols << formatBaselineParams(params, fixed, rms, masklist, whichrow, false) << LogIO::POST ; } if (outTextFile) { ofs << formatBaselineParams(params, fixed, rms, masklist, whichrow, true) << flush; } } } /* for cspline. will be merged once cspline is available in fitter (2011/3/10 WK) */ void Scantable::outputFittingResult(bool outLogger, bool outTextFile, const std::vector& chanMask, int whichrow, const casa::String& coordInfo, bool hasSameNchan, ofstream& ofs, const casa::String& funcName, const std::vector& edge, const std::vector& params, const std::vector& fixed) { if (outLogger || outTextFile) { float rms = getRms(chanMask, whichrow); String masklist = getMaskRangeList(chanMask, whichrow, coordInfo, hasSameNchan); if (outLogger) { LogIO ols(LogOrigin("Scantable", funcName, WHERE)); ols << formatPiecewiseBaselineParams(edge, params, fixed, rms, masklist, whichrow, false) << LogIO::POST ; } if (outTextFile) { ofs << formatPiecewiseBaselineParams(edge, params, fixed, rms, masklist, whichrow, true) << flush; } } } /* for sinusoid. will be merged once sinusoid is available in fitter (2011/3/10 WK) */ void Scantable::outputFittingResult(bool outLogger, bool outTextFile, const std::vector& chanMask, int whichrow, const casa::String& coordInfo, bool hasSameNchan, ofstream& ofs, const casa::String& funcName, const std::vector& params, const std::vector& fixed) { if (outLogger || outTextFile) { float rms = getRms(chanMask, whichrow); String masklist = getMaskRangeList(chanMask, whichrow, coordInfo, hasSameNchan); if (outLogger) { LogIO ols(LogOrigin("Scantable", funcName, WHERE)); ols << formatBaselineParams(params, fixed, rms, masklist, whichrow, false) << LogIO::POST ; } if (outTextFile) { ofs << formatBaselineParams(params, fixed, rms, masklist, whichrow, true) << flush; } } } float Scantable::getRms(const std::vector& mask, int whichrow) { Vector spec; specCol_.get(whichrow, spec); float mean = 0.0; float smean = 0.0; int n = 0; for (uInt i = 0; i < spec.nelements(); ++i) { if (mask[i]) { mean += spec[i]; smean += spec[i]*spec[i]; n++; } } mean /= (float)n; smean /= (float)n; return sqrt(smean - mean*mean); } std::string Scantable::formatBaselineParamsHeader(int whichrow, const std::string& masklist, bool verbose) const { ostringstream oss; if (verbose) { for (int i = 0; i < 60; ++i) { oss << "-"; } oss << endl; oss << " Scan[" << getScan(whichrow) << "]"; oss << " Beam[" << getBeam(whichrow) << "]"; oss << " IF[" << getIF(whichrow) << "]"; oss << " Pol[" << getPol(whichrow) << "]"; oss << " Cycle[" << getCycle(whichrow) << "]: " << endl; oss << "Fitter range = " << masklist << endl; oss << "Baseline parameters" << endl; oss << flush; } return String(oss); } std::string Scantable::formatBaselineParamsFooter(float rms, bool verbose) const { ostringstream oss; if (verbose) { oss << "Results of baseline fit" << endl; oss << " rms = " << setprecision(6) << rms << endl; } return String(oss); } std::string Scantable::formatBaselineParams(const std::vector& params, const std::vector& fixed, float rms, const std::string& masklist, int whichrow, bool verbose, int start, int count, bool resetparamid) const { int nParam = (int)(params.size()); if (nParam < 1) { return(" Not fitted"); } else { ostringstream oss; oss << formatBaselineParamsHeader(whichrow, masklist, verbose); if (start < 0) start = 0; if (count < 0) count = nParam; int end = start + count; if (end > nParam) end = nParam; int paramidoffset = (resetparamid) ? (-start) : 0; for (int i = start; i < end; ++i) { if (i > start) { oss << ","; } std::string sFix = ((fixed.size() > 0) && (fixed[i]) && verbose) ? "(fixed)" : ""; oss << " p" << (i+paramidoffset) << sFix << "= " << right << setw(13) << setprecision(6) << params[i]; } oss << endl; oss << formatBaselineParamsFooter(rms, verbose); return String(oss); } } std::string Scantable::formatPiecewiseBaselineParams(const std::vector& ranges, const std::vector& params, const std::vector& fixed, float rms, const std::string& masklist, int whichrow, bool verbose) const { int nOutParam = (int)(params.size()); int nPiece = (int)(ranges.size()) - 1; if (nOutParam < 1) { return(" Not fitted"); } else if (nPiece < 0) { return formatBaselineParams(params, fixed, rms, masklist, whichrow, verbose); } else if (nPiece < 1) { return(" Bad count of the piece edge info"); } else if (nOutParam % nPiece != 0) { return(" Bad count of the output baseline parameters"); } else { int nParam = nOutParam / nPiece; ostringstream oss; oss << formatBaselineParamsHeader(whichrow, masklist, verbose); stringstream ss; ss << ranges[nPiece] << flush; int wRange = ss.str().size() * 2 + 5; for (int i = 0; i < nPiece; ++i) { ss.str(""); ss << " [" << ranges[i] << "," << (ranges[i+1]-1) << "]"; oss << left << setw(wRange) << ss.str(); oss << formatBaselineParams(params, fixed, rms, masklist, whichrow, false, i*nParam, nParam, true); } oss << formatBaselineParamsFooter(rms, verbose); return String(oss); } } bool Scantable::hasSameNchanOverIFs() { int nIF = nif(-1); int nCh; int totalPositiveNChan = 0; int nPositiveNChan = 0; for (int i = 0; i < nIF; ++i) { nCh = nchan(i); if (nCh > 0) { totalPositiveNChan += nCh; nPositiveNChan++; } } return (totalPositiveNChan == (nPositiveNChan * nchan(0))); } std::string Scantable::getMaskRangeList(const std::vector& mask, int whichrow, const casa::String& coordInfo, bool hasSameNchan, bool verbose) { if (mask.size() < 2) { throw(AipsError("The mask elements should be > 1")); } int IF = getIF(whichrow); if (mask.size() != (uInt)nchan(IF)) { throw(AipsError("Number of channels in scantable != number of mask elements")); } if (verbose) { LogIO logOs(LogOrigin("Scantable", "getMaskRangeList()", WHERE)); logOs << LogIO::WARN << "The current mask window unit is " << coordInfo; if (!hasSameNchan) { logOs << endl << "This mask is only valid for IF=" << IF; } logOs << LogIO::POST; } std::vector abcissa = getAbcissa(whichrow); std::vector edge = getMaskEdgeIndices(mask); ostringstream oss; oss.setf(ios::fixed); oss << setprecision(1) << "["; for (uInt i = 0; i < edge.size(); i+=2) { if (i > 0) oss << ","; oss << "[" << (float)abcissa[edge[i]] << "," << (float)abcissa[edge[i+1]] << "]"; } oss << "]" << flush; return String(oss); } std::vector Scantable::getMaskEdgeIndices(const std::vector& mask) { if (mask.size() < 2) { throw(AipsError("The mask elements should be > 1")); } std::vector out, startIndices, endIndices; int maskSize = mask.size(); startIndices.clear(); endIndices.clear(); if (mask[0]) { startIndices.push_back(0); } for (int i = 1; i < maskSize; ++i) { if ((!mask[i-1]) && mask[i]) { startIndices.push_back(i); } else if (mask[i-1] && (!mask[i])) { endIndices.push_back(i-1); } } if (mask[maskSize-1]) { endIndices.push_back(maskSize-1); } if (startIndices.size() != endIndices.size()) { throw(AipsError("Inconsistent Mask Size: bad data?")); } for (uInt i = 0; i < startIndices.size(); ++i) { if (startIndices[i] > endIndices[i]) { throw(AipsError("Mask start index > mask end index")); } } out.clear(); for (uInt i = 0; i < startIndices.size(); ++i) { out.push_back(startIndices[i]); out.push_back(endIndices[i]); } return out; } /* STFitEntry Scantable::polyBaseline(const std::vector& mask, int order, int rowno) { Fitter fitter = Fitter(); fitter.setExpression("poly", order); std::vector fmask = getMask(rowno); if (fmask.size() != mask.size()) { throw(AipsError("different mask sizes")); } for (int i = 0; i < fmask.size(); ++i) { fmask[i] = fmask[i] && mask[i]; } fitBaseline(fmask, rowno, fitter); setSpectrum(fitter.getResidual(), rowno); return fitter.getFitEntry(); } */ } //namespace asap