source: trunk/src/SDMath.cc@ 231

Last change on this file since 231 was 230, checked in by kil064, 20 years ago

use new MemomryTable version of readAsciiTable
Restucture correction from table more generically
for reuse with opacity corrections in the future

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 34.1 KB
RevLine 
[2]1//#---------------------------------------------------------------------------
2//# SDMath.cc: A collection of single dish mathematical operations
3//#---------------------------------------------------------------------------
4//# Copyright (C) 2004
[125]5//# ATNF
[2]6//#
7//# This program is free software; you can redistribute it and/or modify it
8//# under the terms of the GNU General Public License as published by the Free
9//# Software Foundation; either version 2 of the License, or (at your option)
10//# any later version.
11//#
12//# This program is distributed in the hope that it will be useful, but
13//# WITHOUT ANY WARRANTY; without even the implied warranty of
14//# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General
15//# Public License for more details.
16//#
17//# You should have received a copy of the GNU General Public License along
18//# with this program; if not, write to the Free Software Foundation, Inc.,
19//# 675 Massachusetts Ave, Cambridge, MA 02139, USA.
20//#
21//# Correspondence concerning this software should be addressed as follows:
22//# Internet email: Malte.Marquarding@csiro.au
23//# Postal address: Malte Marquarding,
24//# Australia Telescope National Facility,
25//# P.O. Box 76,
26//# Epping, NSW, 2121,
27//# AUSTRALIA
28//#
29//# $Id:
30//#---------------------------------------------------------------------------
[38]31#include <vector>
32
[81]33#include <casa/aips.h>
34#include <casa/BasicSL/String.h>
35#include <casa/Arrays/IPosition.h>
36#include <casa/Arrays/Array.h>
[130]37#include <casa/Arrays/ArrayIter.h>
38#include <casa/Arrays/VectorIter.h>
[81]39#include <casa/Arrays/ArrayMath.h>
40#include <casa/Arrays/ArrayLogical.h>
41#include <casa/Arrays/MaskedArray.h>
42#include <casa/Arrays/MaskArrMath.h>
43#include <casa/Arrays/MaskArrLogi.h>
[221]44#include <casa/Containers/Block.h>
45#include <casa/Quanta/QC.h>
[177]46#include <casa/Utilities/Assert.h>
[130]47#include <casa/Exceptions.h>
[2]48
[177]49#include <scimath/Mathematics/VectorKernel.h>
50#include <scimath/Mathematics/Convolver.h>
[227]51#include <scimath/Mathematics/InterpolateArray1D.h>
[177]52
[81]53#include <tables/Tables/Table.h>
54#include <tables/Tables/ScalarColumn.h>
55#include <tables/Tables/ArrayColumn.h>
[227]56#include <tables/Tables/ReadAsciiTable.h>
[2]57
[130]58#include <lattices/Lattices/LatticeUtilities.h>
59#include <lattices/Lattices/RebinLattice.h>
[81]60#include <coordinates/Coordinates/SpectralCoordinate.h>
[130]61#include <coordinates/Coordinates/CoordinateSystem.h>
62#include <coordinates/Coordinates/CoordinateUtil.h>
[221]63#include <coordinates/Coordinates/VelocityAligner.h>
[38]64
65#include "MathUtils.h"
[209]66#include "Definitions.h"
[2]67#include "SDContainer.h"
68#include "SDMemTable.h"
69
70#include "SDMath.h"
71
[125]72using namespace casa;
[83]73using namespace asap;
[2]74
[170]75
76SDMath::SDMath()
77{;}
78
[185]79SDMath::SDMath(const SDMath& other)
[170]80{
81
82// No state
83
84}
85
86SDMath& SDMath::operator=(const SDMath& other)
87{
88 if (this != &other) {
89// No state
90 }
91 return *this;
92}
93
[183]94SDMath::~SDMath()
95{;}
[170]96
[183]97
[185]98CountedPtr<SDMemTable> SDMath::average(const Block<CountedPtr<SDMemTable> >& in,
99 const Vector<Bool>& mask, Bool scanAv,
[227]100 const String& weightStr)
[221]101//Bool alignVelocity)
[130]102//
[144]103// Weighted averaging of spectra from one or more Tables.
[130]104//
105{
[221]106 Bool alignVelocity = False;
[2]107
[163]108// Convert weight type
109
110 WeightType wtType = NONE;
[185]111 convertWeightString(wtType, weightStr);
[163]112
[144]113// Create output Table by cloning from the first table
[2]114
[144]115 SDMemTable* pTabOut = new SDMemTable(*in[0],True);
[130]116
[144]117// Setup
[130]118
[144]119 IPosition shp = in[0]->rowAsMaskedArray(0).shape(); // Must not change
120 Array<Float> arr(shp);
121 Array<Bool> barr(shp);
[221]122 const Bool useMask = (mask.nelements() == shp(asap::ChanAxis));
[130]123
[144]124// Columns from Tables
[130]125
[144]126 ROArrayColumn<Float> tSysCol;
127 ROScalarColumn<Double> mjdCol;
128 ROScalarColumn<String> srcNameCol;
129 ROScalarColumn<Double> intCol;
130 ROArrayColumn<uInt> fqIDCol;
[130]131
[144]132// Create accumulation MaskedArray. We accumulate for each channel,if,pol,beam
133// Note that the mask of the accumulation array will ALWAYS remain ALL True.
134// The MA is only used so that when data which is masked Bad is added to it,
135// that data does not contribute.
136
137 Array<Float> zero(shp);
138 zero=0.0;
139 Array<Bool> good(shp);
140 good = True;
141 MaskedArray<Float> sum(zero,good);
142
143// Counter arrays
144
145 Array<Float> nPts(shp); // Number of points
146 nPts = 0.0;
147 Array<Float> nInc(shp); // Increment
148 nInc = 1.0;
149
150// Create accumulation Array for variance. We accumulate for
151// each if,pol,beam, but average over channel. So we need
152// a shape with one less axis dropping channels.
153
154 const uInt nAxesSub = shp.nelements() - 1;
155 IPosition shp2(nAxesSub);
156 for (uInt i=0,j=0; i<(nAxesSub+1); i++) {
[221]157 if (i!=asap::ChanAxis) {
[144]158 shp2(j) = shp(i);
159 j++;
160 }
[2]161 }
[144]162 Array<Float> sumSq(shp2);
163 sumSq = 0.0;
164 IPosition pos2(nAxesSub,0); // For indexing
[130]165
[144]166// Time-related accumulators
[130]167
[144]168 Double time;
169 Double timeSum = 0.0;
170 Double intSum = 0.0;
171 Double interval = 0.0;
[130]172
[144]173// To get the right shape for the Tsys accumulator we need to
174// access a column from the first table. The shape of this
175// array must not change
[130]176
[144]177 Array<Float> tSysSum;
178 {
179 const Table& tabIn = in[0]->table();
180 tSysCol.attach(tabIn,"TSYS");
181 tSysSum.resize(tSysCol.shape(0));
182 }
183 tSysSum =0.0;
184 Array<Float> tSys;
185
186// Scan and row tracking
187
188 Int oldScanID = 0;
189 Int outScanID = 0;
190 Int scanID = 0;
191 Int rowStart = 0;
192 Int nAccum = 0;
193 Int tableStart = 0;
194
195// Source and FreqID
196
197 String sourceName, oldSourceName, sourceNameStart;
198 Vector<uInt> freqID, freqIDStart, oldFreqID;
199
[221]200// Velocity Aligner. We need an aligner for each Direction and FreqID
201// combination. I don't think there is anyway to know how many
202// directions there are.
203// For now, assume all Tables have the same Frequency Table
204
205/*
206 {
207 MEpoch::Ref timeRef(MEpoch::UTC); // Should be in header
208 MDirection::Types dirRef(MDirection::J2000); // Should be in header
209//
210 SDHeader sh = in[0].getSDHeader();
211 const uInt nChan = sh.nchan;
212//
213 const SDFrequencyTable freqTab = in[0]->getSDFreqTable();
214 const uInt nFreqID = freqTab.length();
215 PtrBlock<const VelocityAligner<Float>* > vA(nFreqID);
216
217// Get first time from first table
218
219 const Table& tabIn0 = in[0]->table();
220 mjdCol.attach(tabIn0, "TIME");
221 Double dTmp;
222 mjdCol.get(0, dTmp);
223 MVEpoch tmp2(Quantum<Double>(dTmp, Unit(String("d"))));
224 MEpoch epoch(tmp2, timeRef);
225//
226 for (uInt freqID=0; freqID<nFreqID; freqID++) {
227 SpectralCoordinate sC = in[0]->getCoordinate(freqID);
228 vA[freqID] = new VelocityAligner<Float>(sC, nChan, epoch, const MDirection& dir,
229 const MPosition& pos, const String& velUnit,
230 MDoppler::Types velType, MFrequency::Types velFreqSystem)
231 }
232 }
233*/
234
[144]235// Loop over tables
236
237 Float fac = 1.0;
238 const uInt nTables = in.nelements();
239 for (uInt iTab=0; iTab<nTables; iTab++) {
240
[221]241// Should check that the frequency tables don't change if doing VelocityAlignment
242
[144]243// Attach columns to Table
244
245 const Table& tabIn = in[iTab]->table();
246 tSysCol.attach(tabIn, "TSYS");
247 mjdCol.attach(tabIn, "TIME");
248 srcNameCol.attach(tabIn, "SRCNAME");
249 intCol.attach(tabIn, "INTERVAL");
250 fqIDCol.attach(tabIn, "FREQID");
251
252// Loop over rows in Table
253
254 const uInt nRows = in[iTab]->nRow();
255 for (uInt iRow=0; iRow<nRows; iRow++) {
256
257// Check conformance
258
259 IPosition shp2 = in[iTab]->rowAsMaskedArray(iRow).shape();
260 if (!shp.isEqual(shp2)) {
261 throw (AipsError("Shapes for all rows must be the same"));
262 }
263
264// If we are not doing scan averages, make checks for source and
265// frequency setup and warn if averaging across them
266
267// Get copy of Scan Container for this row
268
269 SDContainer sc = in[iTab]->getSDContainer(iRow);
270 scanID = sc.scanid;
271
272// Get quantities from columns
273
274 srcNameCol.getScalar(iRow, sourceName);
275 mjdCol.get(iRow, time);
276 tSysCol.get(iRow, tSys);
277 intCol.get(iRow, interval);
278 fqIDCol.get(iRow, freqID);
279
280// Initialize first source and freqID
281
282 if (iRow==0 && iTab==0) {
283 sourceNameStart = sourceName;
284 freqIDStart = freqID;
285 }
286
287// If we are doing scan averages, see if we are at the end of an
288// accumulation period (scan). We must check soutce names too,
289// since we might have two tables with one scan each but different
290// source names; we shouldn't average different sources together
291
292 if (scanAv && ( (scanID != oldScanID) ||
293 (iRow==0 && iTab>0 && sourceName!=oldSourceName))) {
294
295// Normalize data in 'sum' accumulation array according to weighting scheme
296
[221]297 normalize(sum, sumSq, nPts, wtType, asap::ChanAxis, nAxesSub);
[144]298
299// Fill scan container. The source and freqID come from the
300// first row of the first table that went into this average (
301// should be the same for all rows in the scan average)
302
303 Float nR(nAccum);
[185]304 fillSDC(sc, sum.getMask(), sum.getArray(), tSysSum/nR, outScanID,
[144]305 timeSum/nR, intSum, sourceNameStart, freqIDStart);
306
307// Write container out to Table
308
309 pTabOut->putSDContainer(sc);
310
311// Reset accumulators
312
313 sum = 0.0;
314 sumSq = 0.0;
315 nAccum = 0;
316//
317 tSysSum =0.0;
318 timeSum = 0.0;
319 intSum = 0.0;
[221]320 nPts = 0.0;
[144]321
322// Increment
323
324 rowStart = iRow; // First row for next accumulation
325 tableStart = iTab; // First table for next accumulation
326 sourceNameStart = sourceName; // First source name for next accumulation
327 freqIDStart = freqID; // First FreqID for next accumulation
328//
329 oldScanID = scanID;
330 outScanID += 1; // Scan ID for next accumulation period
[227]331 }
[144]332
[146]333// Accumulate
[144]334
[185]335 accumulate(timeSum, intSum, nAccum, sum, sumSq, nPts, tSysSum,
[221]336 tSys, nInc, mask, time, interval, in, iTab, iRow, asap::ChanAxis,
[146]337 nAxesSub, useMask, wtType);
[144]338//
339 oldSourceName = sourceName;
340 oldFreqID = freqID;
[184]341 }
[144]342 }
343
344// OK at this point we have accumulation data which is either
345// - accumulated from all tables into one row
346// or
347// - accumulated from the last scan average
348//
349// Normalize data in 'sum' accumulation array according to weighting scheme
[221]350 normalize(sum, sumSq, nPts, wtType, asap::ChanAxis, nAxesSub);
[144]351
352// Create and fill container. The container we clone will be from
353// the last Table and the first row that went into the current
354// accumulation. It probably doesn't matter that much really...
355
356 Float nR(nAccum);
357 SDContainer sc = in[tableStart]->getSDContainer(rowStart);
[185]358 fillSDC(sc, sum.getMask(), sum.getArray(), tSysSum/nR, outScanID,
[144]359 timeSum/nR, intSum, sourceNameStart, freqIDStart);
[221]360 pTabOut->putSDContainer(sc);
[144]361//
362 return CountedPtr<SDMemTable>(pTabOut);
[2]363}
[9]364
[144]365
366
[85]367CountedPtr<SDMemTable>
368SDMath::quotient(const CountedPtr<SDMemTable>& on,
[185]369 const CountedPtr<SDMemTable>& off)
370{
[130]371//
372// Compute quotient spectrum
373//
374 const uInt nRows = on->nRow();
375 if (off->nRow() != nRows) {
376 throw (AipsError("Input Scan Tables must have the same number of rows"));
377 }
[85]378
[130]379// Input Tables and columns
380
[9]381 Table ton = on->table();
382 Table toff = off->table();
[85]383 ROArrayColumn<Float> tsys(toff, "TSYS");
[9]384 ROScalarColumn<Double> mjd(ton, "TIME");
[15]385 ROScalarColumn<Double> integr(ton, "INTERVAL");
[9]386 ROScalarColumn<String> srcn(ton, "SRCNAME");
[38]387 ROArrayColumn<uInt> freqidc(ton, "FREQID");
388
[130]389// Output Table cloned from input
[85]390
[171]391 SDMemTable* pTabOut = new SDMemTable(*on, True);
[130]392
393// Loop over rows
394
395 for (uInt i=0; i<nRows; i++) {
396 MaskedArray<Float> mon(on->rowAsMaskedArray(i));
397 MaskedArray<Float> moff(off->rowAsMaskedArray(i));
398 IPosition ipon = mon.shape();
399 IPosition ipoff = moff.shape();
400//
401 Array<Float> tsarr;
402 tsys.get(i, tsarr);
403 if (ipon != ipoff && ipon != tsarr.shape()) {
404 throw(AipsError("on/off not conformant"));
405 }
406
407// Compute quotient
408
409 MaskedArray<Float> tmp = (mon-moff);
410 Array<Float> out(tmp.getArray());
411 out /= moff;
412 out *= tsarr;
[163]413 Array<Bool> outflagsb = mon.getMask() && moff.getMask();
[130]414
415// Fill container for this row
416
[184]417 SDContainer sc = on->getSDContainer(i);
[163]418//
[185]419 putDataInSDC(sc, out, outflagsb);
[130]420 sc.putTsys(tsarr);
[184]421 sc.scanid = i;
[130]422
423// Put new row in output Table
424
[171]425 pTabOut->putSDContainer(sc);
[130]426 }
427//
[171]428 return CountedPtr<SDMemTable>(pTabOut);
[9]429}
[48]430
[146]431
432
[185]433std::vector<float> SDMath::statistic(const CountedPtr<SDMemTable>& in,
434 const std::vector<bool>& mask,
435 const String& which)
[130]436//
437// Perhaps iteration over pol/beam/if should be in here
438// and inside the nrow iteration ?
439//
440{
441 const uInt nRow = in->nRow();
442 std::vector<float> result(nRow);
443 Vector<Bool> msk(mask);
444
445// Specify cursor location
446
[152]447 IPosition start, end;
[185]448 getCursorLocation(start, end, *in);
[130]449
450// Loop over rows
451
452 const uInt nEl = msk.nelements();
453 for (uInt ii=0; ii < in->nRow(); ++ii) {
454
455// Get row and deconstruct
456
457 MaskedArray<Float> marr(in->rowAsMaskedArray(ii));
458 Array<Float> arr = marr.getArray();
459 Array<Bool> barr = marr.getMask();
460
461// Access desired piece of data
462
463 Array<Float> v((arr(start,end)).nonDegenerate());
464 Array<Bool> m((barr(start,end)).nonDegenerate());
465
466// Apply OTF mask
467
468 MaskedArray<Float> tmp;
469 if (m.nelements()==nEl) {
470 tmp.setData(v,m&&msk);
471 } else {
472 tmp.setData(v,m);
473 }
474
475// Get statistic
476
[144]477 result[ii] = mathutil::statistics(which, tmp);
[130]478 }
479//
480 return result;
481}
482
[146]483
[185]484SDMemTable* SDMath::bin(const SDMemTable& in, Int width)
[144]485{
[169]486 SDHeader sh = in.getSDHeader();
487 SDMemTable* pTabOut = new SDMemTable(in, True);
[163]488
[169]489// Bin up SpectralCoordinates
[163]490
[169]491 IPosition factors(1);
492 factors(0) = width;
493 for (uInt j=0; j<in.nCoordinates(); ++j) {
494 CoordinateSystem cSys;
495 cSys.addCoordinate(in.getCoordinate(j));
496 CoordinateSystem cSysBin =
[185]497 CoordinateUtil::makeBinnedCoordinateSystem(factors, cSys, False);
[169]498//
499 SpectralCoordinate sCBin = cSysBin.spectralCoordinate(0);
500 pTabOut->setCoordinate(sCBin, j);
501 }
[163]502
[169]503// Use RebinLattice to find shape
[130]504
[169]505 IPosition shapeIn(1,sh.nchan);
[185]506 IPosition shapeOut = RebinLattice<Float>::rebinShape(shapeIn, factors);
[169]507 sh.nchan = shapeOut(0);
508 pTabOut->putSDHeader(sh);
[144]509
510
[169]511// Loop over rows and bin along channel axis
512
513 for (uInt i=0; i < in.nRow(); ++i) {
514 SDContainer sc = in.getSDContainer(i);
[144]515//
[169]516 Array<Float> tSys(sc.getTsys()); // Get it out before sc changes shape
[144]517
[169]518// Bin up spectrum
[144]519
[169]520 MaskedArray<Float> marr(in.rowAsMaskedArray(i));
521 MaskedArray<Float> marrout;
[221]522 LatticeUtilities::bin(marrout, marr, asap::ChanAxis, width);
[144]523
[169]524// Put back the binned data and flags
[144]525
[169]526 IPosition ip2 = marrout.shape();
527 sc.resize(ip2);
[146]528//
[185]529 putDataInSDC(sc, marrout.getArray(), marrout.getMask());
[146]530
[169]531// Bin up Tsys.
[146]532
[169]533 Array<Bool> allGood(tSys.shape(),True);
534 MaskedArray<Float> tSysIn(tSys, allGood, True);
[146]535//
[169]536 MaskedArray<Float> tSysOut;
[221]537 LatticeUtilities::bin(tSysOut, tSysIn, asap::ChanAxis, width);
[169]538 sc.putTsys(tSysOut.getArray());
[146]539//
[169]540 pTabOut->putSDContainer(sc);
541 }
542 return pTabOut;
[146]543}
544
[185]545SDMemTable* SDMath::simpleOperate(const SDMemTable& in, Float val, Bool doAll,
546 uInt what)
[152]547//
548// what = 0 Multiply
549// 1 Add
[146]550{
[152]551 SDMemTable* pOut = new SDMemTable(in,False);
552 const Table& tOut = pOut->table();
553 ArrayColumn<Float> spec(tOut,"SPECTRA");
[146]554//
[152]555 if (doAll) {
556 for (uInt i=0; i < tOut.nrow(); i++) {
557
558// Get
559
560 MaskedArray<Float> marr(pOut->rowAsMaskedArray(i));
561
562// Operate
563
564 if (what==0) {
565 marr *= val;
566 } else if (what==1) {
567 marr += val;
568 }
569
570// Put
571
572 spec.put(i, marr.getArray());
573 }
574 } else {
575
576// Get cursor location
577
578 IPosition start, end;
[185]579 getCursorLocation(start, end, in);
[152]580//
581 for (uInt i=0; i < tOut.nrow(); i++) {
582
583// Get
584
585 MaskedArray<Float> dataIn(pOut->rowAsMaskedArray(i));
586
587// Modify. More work than we would like to deal with the mask
588
589 Array<Float>& values = dataIn.getRWArray();
590 Array<Bool> mask(dataIn.getMask());
591//
592 Array<Float> values2 = values(start,end);
593 Array<Bool> mask2 = mask(start,end);
594 MaskedArray<Float> t(values2,mask2);
595 if (what==0) {
596 t *= val;
597 } else if (what==1) {
598 t += val;
599 }
600 values(start, end) = t.getArray(); // Write back into 'dataIn'
601
602// Put
603 spec.put(i, dataIn.getArray());
604 }
605 }
606//
[146]607 return pOut;
608}
609
610
[152]611
[185]612SDMemTable* SDMath::averagePol(const SDMemTable& in, const Vector<Bool>& mask)
[152]613//
[165]614// Average all polarizations together, weighted by variance
615//
616{
617// WeightType wtType = NONE;
[185]618// convertWeightString(wtType, weight);
[165]619
620 const uInt nRows = in.nRow();
[209]621 const uInt polAxis = asap::PolAxis; // Polarization axis
622 const uInt chanAxis = asap::ChanAxis; // Spectrum axis
[165]623
624// Create output Table and reshape number of polarizations
625
626 Bool clear=True;
627 SDMemTable* pTabOut = new SDMemTable(in, clear);
628 SDHeader header = pTabOut->getSDHeader();
629 header.npol = 1;
630 pTabOut->putSDHeader(header);
631
632// Shape of input and output data
633
634 const IPosition& shapeIn = in.rowAsMaskedArray(0u, False).shape();
635 IPosition shapeOut(shapeIn);
636 shapeOut(polAxis) = 1; // Average all polarizations
637//
638 const uInt nChan = shapeIn(chanAxis);
639 const IPosition vecShapeOut(4,1,1,1,nChan); // A multi-dim form of a Vector shape
640 IPosition start(4), end(4);
641
642// Output arrays
643
644 Array<Float> outData(shapeOut, 0.0);
645 Array<Bool> outMask(shapeOut, True);
646 const IPosition axes(2, 2, 3); // pol-channel plane
647//
648 const Bool useMask = (mask.nelements() == shapeIn(chanAxis));
649
650// Loop over rows
651
652 for (uInt iRow=0; iRow<nRows; iRow++) {
653
654// Get data for this row
655
656 MaskedArray<Float> marr(in.rowAsMaskedArray(iRow));
657 Array<Float>& arr = marr.getRWArray();
658 const Array<Bool>& barr = marr.getMask();
659
660// Make iterators to iterate by pol-channel planes
661
662 ReadOnlyArrayIterator<Float> itDataPlane(arr, axes);
663 ReadOnlyArrayIterator<Bool> itMaskPlane(barr, axes);
664
665// Accumulations
666
667 Float fac = 1.0;
668 Vector<Float> vecSum(nChan,0.0);
669
670// Iterate through data by pol-channel planes
671
672 while (!itDataPlane.pastEnd()) {
673
674// Iterate through plane by polarization and accumulate Vectors
675
676 Vector<Float> t1(nChan); t1 = 0.0;
677 Vector<Bool> t2(nChan); t2 = True;
678 MaskedArray<Float> vecSum(t1,t2);
679 Float varSum = 0.0;
680 {
681 ReadOnlyVectorIterator<Float> itDataVec(itDataPlane.array(), 1);
682 ReadOnlyVectorIterator<Bool> itMaskVec(itMaskPlane.array(), 1);
683 while (!itDataVec.pastEnd()) {
684
685// Create MA of data & mask (optionally including OTF mask) and get variance
686
687 if (useMask) {
688 const MaskedArray<Float> spec(itDataVec.vector(),mask&&itMaskVec.vector());
689 fac = 1.0 / variance(spec);
690 } else {
691 const MaskedArray<Float> spec(itDataVec.vector(),itMaskVec.vector());
692 fac = 1.0 / variance(spec);
693 }
694
695// Normalize spectrum (without OTF mask) and accumulate
696
697 const MaskedArray<Float> spec(fac*itDataVec.vector(), itMaskVec.vector());
698 vecSum += spec;
699 varSum += fac;
700
701// Next
702
703 itDataVec.next();
704 itMaskVec.next();
705 }
706 }
707
708// Normalize summed spectrum
709
710 vecSum /= varSum;
711
712// FInd position in input data array. We are iterating by pol-channel
713// plane so all that will change is beam and IF and that's what we want.
714
715 IPosition pos = itDataPlane.pos();
716
717// Write out data. This is a bit messy. We have to reform the Vector
718// accumulator into an Array of shape (1,1,1,nChan)
719
720 start = pos;
721 end = pos;
722 end(chanAxis) = nChan-1;
723 outData(start,end) = vecSum.getArray().reform(vecShapeOut);
724 outMask(start,end) = vecSum.getMask().reform(vecShapeOut);
725
726// Step to next beam/IF combination
727
728 itDataPlane.next();
729 itMaskPlane.next();
730 }
731
732// Generate output container and write it to output table
733
734 SDContainer sc = in.getSDContainer();
735 sc.resize(shapeOut);
736//
[185]737 putDataInSDC(sc, outData, outMask);
[165]738 pTabOut->putSDContainer(sc);
739 }
740//
741 return pTabOut;
742}
[167]743
[169]744
[185]745SDMemTable* SDMath::smooth(const SDMemTable& in,
746 const casa::String& kernelType,
747 casa::Float width, Bool doAll)
[177]748{
[169]749
[177]750// Number of channels
[169]751
[209]752 const uInt chanAxis = asap::ChanAxis; // Spectral axis
[177]753 SDHeader sh = in.getSDHeader();
754 const uInt nChan = sh.nchan;
755
756// Generate Kernel
757
[185]758 VectorKernel::KernelTypes type = VectorKernel::toKernelType(kernelType);
[177]759 Vector<Float> kernel = VectorKernel::make(type, width, nChan, True, False);
760
761// Generate Convolver
762
763 IPosition shape(1,nChan);
764 Convolver<Float> conv(kernel, shape);
765
766// New Table
767
768 SDMemTable* pTabOut = new SDMemTable(in,True);
769
770// Get cursor location
771
772 IPosition start, end;
[185]773 getCursorLocation(start, end, in);
[177]774//
775 IPosition shapeOut(4,1);
776
777// Output Vectors
778
779 Vector<Float> valuesOut(nChan);
780 Vector<Bool> maskOut(nChan);
781
782// Loop over rows in Table
783
784 for (uInt ri=0; ri < in.nRow(); ++ri) {
785
786// Get copy of data
787
788 const MaskedArray<Float>& dataIn(in.rowAsMaskedArray(ri));
789 AlwaysAssert(dataIn.shape()(chanAxis)==nChan, AipsError);
790//
791 Array<Float> valuesIn = dataIn.getArray();
792 Array<Bool> maskIn = dataIn.getMask();
793
794// Branch depending on whether we smooth all locations or just
795// those pointed at by the current selection cursor
796
797 if (doAll) {
[221]798 uInt axis = asap::ChanAxis;
[177]799 VectorIterator<Float> itValues(valuesIn, axis);
800 VectorIterator<Bool> itMask(maskIn, axis);
801 while (!itValues.pastEnd()) {
802
803// Smooth
804 if (kernelType==VectorKernel::HANNING) {
805 mathutil::hanning(valuesOut, maskOut, itValues.vector(), itMask.vector());
806 itMask.vector() = maskOut;
807 } else {
808 mathutil::replaceMaskByZero(itValues.vector(), itMask.vector());
809 conv.linearConv(valuesOut, itValues.vector());
810 }
811//
812 itValues.vector() = valuesOut;
813//
814 itValues.next();
815 itMask.next();
816 }
817 } else {
818
819// Set multi-dim Vector shape
820
821 shapeOut(chanAxis) = valuesIn.shape()(chanAxis);
822
823// Stuff about with shapes so that we don't have conformance run-time errors
824
825 Vector<Float> valuesIn2 = valuesIn(start,end).nonDegenerate();
826 Vector<Bool> maskIn2 = maskIn(start,end).nonDegenerate();
827
828// Smooth
829
830 if (kernelType==VectorKernel::HANNING) {
831 mathutil::hanning(valuesOut, maskOut, valuesIn2, maskIn2);
832 maskIn(start,end) = maskOut.reform(shapeOut);
833 } else {
834 mathutil::replaceMaskByZero(valuesIn2, maskIn2);
835 conv.linearConv(valuesOut, valuesIn2);
836 }
837//
838 valuesIn(start,end) = valuesOut.reform(shapeOut);
839 }
840
841// Create and put back
842
843 SDContainer sc = in.getSDContainer(ri);
[185]844 putDataInSDC(sc, valuesIn, maskIn);
[177]845//
846 pTabOut->putSDContainer(sc);
847 }
848//
849 return pTabOut;
850}
851
852
[221]853SDMemTable* SDMath::convertFlux (const SDMemTable& in, Float a, Float eta, Bool doAll)
854//
855// As it is, this function could be implemented with 'simpleOperate'
856// However, I anticipate that eventually we will look the conversion
857// values up in a Table and apply them in a frequency dependent way,
858// so I have implemented it fully here
859//
860{
861 SDHeader sh = in.getSDHeader();
862 SDMemTable* pTabOut = new SDMemTable(in, True);
[177]863
[221]864// FInd out how to convert values into Jy and K (e.g. units might be mJy or mK)
865// Also automatically find out what we are converting to according to the
866// flux unit
[177]867
[221]868 Unit fluxUnit(sh.fluxunit);
869 Unit K(String("K"));
870 Unit JY(String("Jy"));
871//
872 Bool toKelvin = True;
873 Double inFac = 1.0;
874 if (fluxUnit==JY) {
875 cerr << "Converting to K" << endl;
876//
877 Quantum<Double> t(1.0,fluxUnit);
878 Quantum<Double> t2 = t.get(JY);
879 inFac = (t2 / t).getValue();
880//
881 toKelvin = True;
882 sh.fluxunit = "K";
883 } else if (fluxUnit==K) {
884 cerr << "Converting to Jy" << endl;
885//
886 Quantum<Double> t(1.0,fluxUnit);
887 Quantum<Double> t2 = t.get(K);
888 inFac = (t2 / t).getValue();
889//
890 toKelvin = False;
891 sh.fluxunit = "Jy";
892 } else {
893 throw AipsError("Unrecognized brightness units in Table - must be consistent with Jy or K");
894 }
895 pTabOut->putSDHeader(sh);
[177]896
[221]897// Compute conversion factor. 'a' and 'eta' are really frequency, time and
898// telescope dependent and should be looked// up in a table
899
900 Float factor = 2.0 * inFac * 1.0e-7 * 1.0e26 * QC::k.getValue(Unit(String("erg/K"))) / a / eta;
901 if (toKelvin) {
902 factor = 1.0 / factor;
903 }
904 cerr << "Applying conversion factor = " << factor << endl;
905
906// For operations only on specified cursor location
907
908 IPosition start, end;
909 getCursorLocation(start, end, in);
910
911// Loop over rows and apply factor to spectra
912
913 const uInt axis = asap::ChanAxis;
914 for (uInt i=0; i < in.nRow(); ++i) {
915
916// Get data
917
918 MaskedArray<Float> dataIn(in.rowAsMaskedArray(i));
919 Array<Float>& valuesIn = dataIn.getRWArray(); // writable reference
920 const Array<Bool>& maskIn = dataIn.getMask();
921
922// Need to apply correct conversion factor (frequency and time dependent)
923// which should be sourced from a Table. For now we just apply the given
924// factor to everything
925
926 if (doAll) {
927 VectorIterator<Float> itValues(valuesIn, asap::ChanAxis);
928 while (!itValues.pastEnd()) {
929 itValues.vector() *= factor; // Writes back into dataIn
930//
931 itValues.next();
932 }
933 } else {
934 Array<Float> valuesIn2 = valuesIn(start,end);
935 valuesIn2 *= factor;
936 valuesIn(start,end) = valuesIn2;
937 }
938
939// Write out
940
941 SDContainer sc = in.getSDContainer(i);
942 putDataInSDC(sc, valuesIn, maskIn);
943//
944 pTabOut->putSDContainer(sc);
945 }
946 return pTabOut;
947}
948
949
950
[227]951SDMemTable* SDMath::gainElevation (const SDMemTable& in, const String& fileName,
952 const String& methodStr, Bool doAll)
953{
954 SDHeader sh = in.getSDHeader();
955 SDMemTable* pTabOut = new SDMemTable(in, True);
956 const uInt nRow = in.nRow();
957
[230]958// Get elevation from SDMemTable data
[227]959
960 const Table& tab = in.table();
961 ROScalarColumn<Float> elev(tab, "ELEVATION");
962 Vector<Float> xOut = elev.getColumn();
963 xOut *= Float(180 / C::pi);
964//
[230]965 String col0("ELEVATION");
966 String col1("FACTOR");
[227]967//
[230]968 return correctFromAsciiTable (pTabOut, in, fileName, col0, col1, methodStr, doAll, xOut);
969}
[227]970
[230]971
[227]972
973
974
[169]975// 'private' functions
976
[185]977void SDMath::fillSDC(SDContainer& sc,
978 const Array<Bool>& mask,
979 const Array<Float>& data,
980 const Array<Float>& tSys,
981 Int scanID, Double timeStamp,
982 Double interval, const String& sourceName,
[227]983 const Vector<uInt>& freqID) const
[167]984{
[169]985// Data and mask
[167]986
[185]987 putDataInSDC(sc, data, mask);
[167]988
[169]989// TSys
990
991 sc.putTsys(tSys);
992
993// Time things
994
995 sc.timestamp = timeStamp;
996 sc.interval = interval;
997 sc.scanid = scanID;
[167]998//
[169]999 sc.sourcename = sourceName;
1000 sc.putFreqMap(freqID);
1001}
[167]1002
[185]1003void SDMath::normalize(MaskedArray<Float>& sum,
[169]1004 const Array<Float>& sumSq,
1005 const Array<Float>& nPts,
1006 WeightType wtType, Int axis,
[227]1007 Int nAxesSub) const
[169]1008{
1009 IPosition pos2(nAxesSub,0);
1010//
1011 if (wtType==NONE) {
[167]1012
[169]1013// We just average by the number of points accumulated.
1014// We need to make a MA out of nPts so that no divide by
1015// zeros occur
[167]1016
[169]1017 MaskedArray<Float> t(nPts, (nPts>Float(0.0)));
1018 sum /= t;
1019 } else if (wtType==VAR) {
[167]1020
[169]1021// Normalize each spectrum by sum(1/var) where the variance
1022// is worked out for each spectrum
1023
1024 Array<Float>& data = sum.getRWArray();
1025 VectorIterator<Float> itData(data, axis);
1026 while (!itData.pastEnd()) {
1027 pos2 = itData.pos().getFirst(nAxesSub);
1028 itData.vector() /= sumSq(pos2);
1029 itData.next();
1030 }
1031 } else if (wtType==TSYS) {
1032 }
1033}
1034
1035
[185]1036void SDMath::accumulate(Double& timeSum, Double& intSum, Int& nAccum,
1037 MaskedArray<Float>& sum, Array<Float>& sumSq,
1038 Array<Float>& nPts, Array<Float>& tSysSum,
1039 const Array<Float>& tSys, const Array<Float>& nInc,
1040 const Vector<Bool>& mask, Double time, Double interval,
1041 const Block<CountedPtr<SDMemTable> >& in,
1042 uInt iTab, uInt iRow, uInt axis,
1043 uInt nAxesSub, Bool useMask,
[227]1044 WeightType wtType) const
[169]1045{
1046
1047// Get data
1048
1049 MaskedArray<Float> dataIn(in[iTab]->rowAsMaskedArray(iRow));
1050 Array<Float>& valuesIn = dataIn.getRWArray(); // writable reference
1051 const Array<Bool>& maskIn = dataIn.getMask(); // RO reference
[167]1052//
[169]1053 if (wtType==NONE) {
1054 const MaskedArray<Float> n(nInc,dataIn.getMask());
1055 nPts += n; // Only accumulates where mask==T
1056 } else if (wtType==VAR) {
[167]1057
[169]1058// We are going to average the data, weighted by the noise for each pol, beam and IF.
1059// So therefore we need to iterate through by spectrum (axis 3)
[167]1060
[169]1061 VectorIterator<Float> itData(valuesIn, axis);
1062 ReadOnlyVectorIterator<Bool> itMask(maskIn, axis);
1063 Float fac = 1.0;
1064 IPosition pos(nAxesSub,0);
1065//
1066 while (!itData.pastEnd()) {
[167]1067
[169]1068// Make MaskedArray of Vector, optionally apply OTF mask, and find scaling factor
[167]1069
[169]1070 if (useMask) {
1071 MaskedArray<Float> tmp(itData.vector(),mask&&itMask.vector());
1072 fac = 1.0/variance(tmp);
1073 } else {
1074 MaskedArray<Float> tmp(itData.vector(),itMask.vector());
1075 fac = 1.0/variance(tmp);
1076 }
1077
1078// Scale data
1079
1080 itData.vector() *= fac; // Writes back into 'dataIn'
[167]1081//
[169]1082// Accumulate variance per if/pol/beam averaged over spectrum
1083// This method to get pos2 from itData.pos() is only valid
1084// because the spectral axis is the last one (so we can just
1085// copy the first nAXesSub positions out)
[167]1086
[169]1087 pos = itData.pos().getFirst(nAxesSub);
1088 sumSq(pos) += fac;
1089//
1090 itData.next();
1091 itMask.next();
1092 }
1093 } else if (wtType==TSYS) {
1094 }
[167]1095
[169]1096// Accumulate sum of (possibly scaled) data
1097
1098 sum += dataIn;
1099
1100// Accumulate Tsys, time, and interval
1101
1102 tSysSum += tSys;
1103 timeSum += time;
1104 intSum += interval;
1105 nAccum += 1;
1106}
1107
1108
1109
1110
[185]1111void SDMath::getCursorLocation(IPosition& start, IPosition& end,
[227]1112 const SDMemTable& in) const
[169]1113{
1114 const uInt nDim = 4;
1115 const uInt i = in.getBeam();
1116 const uInt j = in.getIF();
1117 const uInt k = in.getPol();
1118 const uInt n = in.nChan();
[167]1119//
[169]1120 start.resize(nDim);
1121 start(0) = i;
1122 start(1) = j;
1123 start(2) = k;
1124 start(3) = 0;
[167]1125//
[169]1126 end.resize(nDim);
1127 end(0) = i;
1128 end(1) = j;
1129 end(2) = k;
1130 end(3) = n-1;
1131}
1132
1133
[227]1134void SDMath::convertWeightString(WeightType& wtType, const String& weightStr) const
[169]1135{
1136 String tStr(weightStr);
1137 tStr.upcase();
1138 if (tStr.contains(String("NONE"))) {
1139 wtType = NONE;
1140 } else if (tStr.contains(String("VAR"))) {
1141 wtType = VAR;
1142 } else if (tStr.contains(String("TSYS"))) {
1143 wtType = TSYS;
[185]1144 throw(AipsError("T_sys weighting not yet implemented"));
[169]1145 } else {
[185]1146 throw(AipsError("Unrecognized weighting type"));
[167]1147 }
1148}
1149
[227]1150void SDMath::convertInterpString(Int& type, const String& interp) const
1151{
1152 String tStr(interp);
1153 tStr.upcase();
1154 if (tStr.contains(String("NEAR"))) {
1155 type = InterpolateArray1D<Float,Float>::nearestNeighbour;
1156 } else if (tStr.contains(String("LIN"))) {
1157 type = InterpolateArray1D<Float,Float>::linear;
1158 } else if (tStr.contains(String("CUB"))) {
1159 type = InterpolateArray1D<Float,Float>::cubic;
1160 } else if (tStr.contains(String("SPL"))) {
1161 type = InterpolateArray1D<Float,Float>::spline;
1162 } else {
1163 throw(AipsError("Unrecognized interpolation type"));
1164 }
1165}
1166
[185]1167void SDMath::putDataInSDC(SDContainer& sc, const Array<Float>& data,
[227]1168 const Array<Bool>& mask) const
[169]1169{
1170 sc.putSpectrum(data);
1171//
1172 Array<uChar> outflags(data.shape());
1173 convertArray(outflags,!mask);
1174 sc.putFlags(outflags);
1175}
[227]1176
1177Table SDMath::readAsciiFile (const String& fileName) const
1178{
[230]1179 String formatString;
1180 Table tbl = readAsciiTable (formatString, Table::Memory, fileName, "", "", False);
[227]1181 return tbl;
1182}
[230]1183
1184
1185SDMemTable* SDMath::correctFromAsciiTable(SDMemTable* pTabOut,
1186 const SDMemTable& in, const String& fileName,
1187 const String& col0, const String& col1,
1188 const String& methodStr, Bool doAll,
1189 const Vector<Float>& xOut)
1190{
1191
1192// Read gain-elevation ascii file data into a Table.
1193
1194 Table tTable = readAsciiFile (fileName);
1195// tTable.markForDelete();
1196//
1197 return correctFromTable (pTabOut, in, tTable, col0, col1, methodStr, doAll, xOut);
1198}
1199
1200SDMemTable* SDMath::correctFromTable(SDMemTable* pTabOut, const SDMemTable& in, const Table& tTable,
1201 const String& col0, const String& col1,
1202 const String& methodStr, Bool doAll,
1203 const Vector<Float>& xOut)
1204{
1205
1206// Get data from Table
1207
1208 ROScalarColumn<Float> geElCol(tTable, col0);
1209 ROScalarColumn<Float> geFacCol(tTable, col1);
1210 Vector<Float> xIn = geElCol.getColumn();
1211 Vector<Float> yIn = geFacCol.getColumn();
1212 Vector<Bool> maskIn(xIn.nelements(),True);
1213
1214// Interpolate (and extrapolate) with desired method
1215
1216 Int method = 0;
1217 convertInterpString(method, methodStr);
1218//
1219 Vector<Float> yOut;
1220 Vector<Bool> maskOut;
1221 InterpolateArray1D<Float,Float>::interpolate(yOut, maskOut, xOut,
1222 xIn, yIn, maskIn, method,
1223 True, True);
1224
1225// For operations only on specified cursor location
1226
1227 IPosition start, end;
1228 getCursorLocation(start, end, in);
1229
1230// Loop over rows and interpolate correction factor
1231
1232 const uInt axis = asap::ChanAxis;
1233 for (uInt i=0; i < in.nRow(); ++i) {
1234
1235// Get data
1236
1237 MaskedArray<Float> dataIn(in.rowAsMaskedArray(i));
1238 Array<Float>& valuesIn = dataIn.getRWArray(); // writable reference
1239 const Array<Bool>& maskIn = dataIn.getMask();
1240
1241// Apply factor
1242
1243 if (doAll) {
1244 VectorIterator<Float> itValues(valuesIn, asap::ChanAxis);
1245 while (!itValues.pastEnd()) {
1246 itValues.vector() *= yOut(i); // Writes back into dataIn
1247//
1248 itValues.next();
1249 }
1250 } else {
1251 Array<Float> valuesIn2 = valuesIn(start,end);
1252 valuesIn2 *= yOut(i);
1253 valuesIn(start,end) = valuesIn2;
1254 }
1255
1256// Write out
1257
1258 SDContainer sc = in.getSDContainer(i);
1259 putDataInSDC(sc, valuesIn, maskIn);
1260//
1261 pTabOut->putSDContainer(sc);
1262 }
1263//
1264 return pTabOut;
1265}
1266
Note: See TracBrowser for help on using the repository browser.