source: trunk/src/STFITSImageWriter.cpp@ 1452

Last change on this file since 1452 was 1444, checked in by Malte Marquarding, 16 years ago

Ticket #110: I have added a fix for the FITS header now, This modifies CRVAL1, CTYPE1 to be read by CLASS

File size: 7.5 KB
Line 
1//#---------------------------------------------------------------------------
2//# STFITSImageWriter.cc: Class to write out single dish spectra as FITS images
3//#---------------------------------------------------------------------------
4//# Copyright (C) 2008
5//# ATNF
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//#---------------------------------------------------------------------------
31
32#include <fitsio.h>
33#include <images/Images/TempImage.h>
34
35#include <lattices/Lattices/ArrayLattice.h>
36
37#include <measures/Measures/MEpoch.h>
38#include <measures/Measures/Stokes.h>
39
40#include <tables/Tables/Table.h>
41#include <tables/Tables/TableIter.h>
42#include <tables/Tables/TableRecord.h>
43#include <casa/Containers/RecordField.h>
44#include <tables/Tables/TableRow.h>
45#include <tables/Tables/ScalarColumn.h>
46#include <tables/Tables/ArrayColumn.h>
47
48
49#include <coordinates/Coordinates/CoordinateUtil.h>
50#include <coordinates/Coordinates/CoordinateSystem.h>
51#include <coordinates/Coordinates/SpectralCoordinate.h>
52#include <coordinates/Coordinates/DirectionCoordinate.h>
53#include <coordinates/Coordinates/StokesCoordinate.h>
54#include <coordinates/Coordinates/ObsInfo.h>
55
56#include <images/Images/ImageFITSConverter.h>
57#include <images/Images/TempImage.h>
58
59#include <lattices/Lattices/ArrayLattice.h>
60
61
62//#include "STDefs.h"
63#include "STHeader.h"
64#include "Scantable.h"
65#include "STFITSImageWriter.h"
66
67using namespace casa;
68using namespace asap;
69
70
71STFITSImageWriter::STFITSImageWriter() : isClass_(False)
72{;}
73
74STFITSImageWriter::~STFITSImageWriter()
75{;}
76
77
78Bool STFITSImageWriter::write(const Scantable& stable,
79 const String& fileName)
80{
81
82// Get global Header from Table
83
84 STHeader hdr = stable.getHeader();
85
86// Column keywords
87
88 Table tab = stable.table();
89
90// Temps
91
92 const Unit RAD(String("rad"));
93
94// Open and write header file
95
96 String rootName(fileName);
97 if (rootName.length()==0) rootName = String("fits");
98
99 ObsInfo oi;
100 oi.setTelescope(String(hdr.antennaname));
101 oi.setObserver(String(hdr.observer));
102
103 uInt maxMem = 128;
104 Bool preferVelocity = False;
105 Bool opticalVelocity = False;
106 Int bitPix = -32; // FLoating point
107 Float minPix = 1.0;
108 Float maxPix = -1.0;
109 Bool overWrite = True;
110 Bool degLast = False;
111 Bool reallyVerbose = False;
112 String errMsg;
113
114
115 Block<String> cols(5);
116 cols[0] = String("SCANNO");
117 cols[1] = String("CYCLENO");
118 cols[2] = String("BEAMNO");
119 cols[3] = String("IFNO");
120 cols[4] = String("POLNO");
121 TableIterator iter(tab, cols);
122 // Open data file
123 while ( !iter.pastEnd() ) {
124 Table t = iter.table();
125 ROTableRow row(t);
126 const TableRecord& rec = row.get(0);
127 String dirtype = stable.getDirectionRefString();
128 ostringstream onstr;
129 onstr << "SCAN" << rec.asuInt("SCANNO")
130 << "_CYCLE" << rec.asuInt("CYCLENO")
131 << "_BEAM" << rec.asuInt("BEAMNO")
132 << "_IF" << rec.asuInt("IFNO")
133 << "_POL" << rec.asuInt("POLNO");
134 String fileName = rootName + String(onstr) + String(".fits");
135 int row0 = t.rowNumbers(tab)[0];
136
137 const MPosition& mp = stable.getAntennaPosition();
138 const MDirection& md = stable.getDirection(row0);
139 const MEpoch& me = stable.getEpoch(row0);
140 oi.setObsDate(me);
141
142 const Double& rf =
143 stable.molecules().getRestFrequency(rec.asuInt("MOLECULE_ID") );
144 SpectralCoordinate sC =
145 stable.frequencies().getSpectralCoordinate(md, mp, me, rf,
146 rec.asuInt("FREQ_ID"));
147
148 Int polno = rec.asuInt("POLNO");
149 Stokes::StokesTypes stokes =
150 Stokes::type(stable.getPolarizationLabel(polno, stable.getPolType()));
151 Vector<Int> whichStokes(1);
152 whichStokes(0) = Int(stokes);
153 StokesCoordinate stC(whichStokes);
154
155
156 CoordinateSystem cSys;
157 Vector<Double> lonlat = md.getAngle().getValue();
158 DirectionCoordinate dC =
159 getDirectionCoordinate(stable.getDirectionRefString(),
160 lonlat[0], lonlat[1]);
161 cSys.addCoordinate(sC);
162 cSys.addCoordinate(dC);
163 cSys.addCoordinate(stC);
164 cSys.setObsInfo(oi);
165
166 Vector<Float> spec = rec.asArrayFloat("SPECTRA");
167 Vector<uChar> flag = rec.asArrayuChar("FLAGTRA");
168 Vector<Bool> bflag(flag.shape());
169 convertArray(bflag, flag);
170
171 // Create casacore Image
172 IPosition shp(4, 1);
173 shp(0)= spec.nelements();
174 cout << shp << spec << endl;
175 TempImage<Float> tIm(shp, cSys);
176 tIm.put(spec);
177 ArrayLattice<Bool> latMask(shp);
178 IPosition where(4,0);
179 IPosition stride(4,1);
180 latMask.putSlice(!bflag, where, stride);
181 tIm.attachMask(latMask);
182 Bool ok = ImageFITSConverter::ImageToFITS(errMsg, tIm, fileName,
183 maxMem, preferVelocity,
184 opticalVelocity, bitPix,
185 minPix, maxPix, overWrite,
186 degLast, reallyVerbose);
187 if (!ok) {
188 throw(AipsError(errMsg));
189 }
190 if (isClass_) {
191 classHackHeader(fileName);
192 }
193 //pushLog(String(oss));
194 ++iter;
195 }
196 return True;
197}
198
199DirectionCoordinate
200STFITSImageWriter::getDirectionCoordinate(const String& reff,
201 Double lon, Double lat)
202{
203
204 Projection proj(Projection::SIN); // What should we use ?
205 Matrix<Double> xForm(2,2);
206 xForm = 0.0; xForm.diagonal() = 1.0;
207 Vector<Double> incLonLat(2,0.0);
208 Vector<Double> refPixLonLat(2,0.0);
209 MDirection::Types mdt;
210 if (!MDirection::getType(mdt, reff)) {
211 throw(AipsError("Illegal Direction frame."));
212 }
213 return DirectionCoordinate(mdt, proj, lon, lat,
214 incLonLat[0], incLonLat[1], xForm,
215 refPixLonLat[0], refPixLonLat[1]);
216
217}
218
219void STFITSImageWriter::classHackHeader(const String& filename) {
220 int status = 0;
221 fitsfile *fptr;
222 cout << "filename" << endl;
223 if( fits_open_file(&fptr, filename.c_str(), READWRITE, &status) )
224 throw AipsError("FCoudn't open fits file for CLASS modification");
225
226 if ( fits_update_key(fptr, TSTRING, "CTYPE1", (char *)"FREQ",
227 NULL, &status) )
228 throw AipsError("Couldn't modify CTYPE1.");
229 float restf,refval,newfreq;
230 fits_read_key(fptr, TFLOAT, "CRVAL1",
231 &refval, NULL, &status);
232 fits_read_key(fptr, TFLOAT, "RESTFREQ",
233 &restf, NULL, &status);
234 newfreq = refval - restf;
235 if ( fits_update_key(fptr, TFLOAT, "CRVAL1", &newfreq, NULL, &status) )
236 throw AipsError("Couldn't modify CRVAL1");
237 fits_close_file(fptr, &status);
238
239}
Note: See TracBrowser for help on using the repository browser.