Line data Source code
1 : //# Imager.h: Imager functionality sits here;
2 : //# Copyright (C) 1996,1997,1998,1999,2000,2001,2002,2003
3 : //# Associated Universities, Inc. Washington DC, USA.
4 : //#
5 : //# This library is free software; you can redistribute it and/or modify it
6 : //# under the terms of the GNU Library General Public License as published by
7 : //# the Free Software Foundation; either version 2 of the License, or (at your
8 : //# option) any later version.
9 : //#
10 : //# This library is distributed in the hope that it will be useful, but WITHOUT
11 : //# ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
12 : //# FITNESS FOR A PARTICULAR PURPOSE. See the GNU Library General Public
13 : //# License for more details.
14 : //#
15 : //# You should have received a copy of the GNU Library General Public License
16 : //# along with this library; if not, write to the Free Software Foundation,
17 : //# Inc., 675 Massachusetts Ave, Cambridge, MA 02139, USA.
18 : //#
19 : //# Correspondence concerning AIPS++ should be addressed as follows:
20 : //# Internet email: casa-feedback@nrao.edu.
21 : //# Postal address: AIPS++ Project Office
22 : //# National Radio Astronomy Observatory
23 : //# 520 Edgemont Road
24 : //#
25 : //# $Id$
26 :
27 : #ifndef SYNTHESIS_IMAGER_H
28 : #define SYNTHESIS_IMAGER_H
29 :
30 : #include <casacore/casa/aips.h>
31 : #include <casacore/casa/OS/Timer.h>
32 : #include <casacore/casa/Containers/Record.h>
33 : #include <casacore/ms/MeasurementSets/MeasurementSet.h>
34 : #include <casacore/casa/Arrays/IPosition.h>
35 : #include <casacore/casa/Quanta/Quantum.h>
36 : #include <components/ComponentModels/ConstantSpectrum.h>
37 :
38 : #include <casacore/measures/Measures/MDirection.h>
39 : #include <components/ComponentModels/FluxStandard.h>
40 : //#include <measures/Measures/MDirection.h
41 : #include <casacore/measures/Measures/MPosition.h>
42 : #include <casacore/measures/Measures/MRadialVelocity.h>
43 :
44 : #include <synthesis/TransformMachines/FTMachine.h>
45 : #include <synthesis/TransformMachines/StokesImageUtil.h>
46 :
47 : #include <synthesis/MeasurementComponents/CleanImageSkyModel.h>
48 : #include <synthesis/TransformMachines/EVLAAperture.h>
49 : #include <synthesis/TransformMachines/BeamSquint.h>
50 : #include <synthesis/MeasurementComponents/WFCleanImageSkyModel.h>
51 : #include <synthesis/MeasurementComponents/ClarkCleanImageSkyModel.h>
52 : #include <synthesis/MeasurementEquations/SkyEquation.h>
53 : #include <synthesis/TransformMachines/ATerm.h>
54 : #include <graphics/GenericPlotter/SimplePlotter.h>
55 :
56 :
57 : namespace casacore{
58 :
59 : class MSHistoryHandler;
60 : class MeasurementSet;
61 : class MDirection;
62 : class MFrequency;
63 : class File;
64 : template<class T> class ImageInterface;
65 : }
66 :
67 : namespace casa { //# NAMESPACE CASA - BEGIN
68 :
69 : // Forward declarations
70 : class VisSet;
71 : class VisImagingWeight_p;
72 : class PBMath;
73 : class VPSkyJones;
74 : class EPJones;
75 : // <summary> Class that contains functions needed for imager </summary>
76 :
77 :
78 : class Imager
79 : {
80 : public:
81 : // Default constructor
82 :
83 : Imager();
84 :
85 : Imager(casacore::MeasurementSet& ms, casacore::Bool compress=false, casacore::Bool useModel=false);
86 : Imager(casacore::MeasurementSet& ms, casacore::Bool compress=false);
87 :
88 : // Copy constructor and assignment operator
89 : Imager(const Imager&);
90 : Imager& operator=(const Imager&);
91 :
92 : // Destructor
93 : virtual ~Imager();
94 :
95 : // open all the subtables as userNoReadLock
96 : virtual casacore::Bool openSubTables();
97 :
98 :
99 : // Lock the ms and its subtables
100 : virtual casacore::Bool lock();
101 :
102 : // Unlock the ms and its subtables
103 : virtual casacore::Bool unlock();
104 :
105 :
106 : // Utility function to do channel selection
107 :
108 : casacore::Bool selectDataChannel(casacore::Vector<casacore::Int>& spectralwindowids,
109 : casacore::String& dataMode,
110 : casacore::Vector<casacore::Int>& dataNchan,
111 : casacore::Vector<casacore::Int>& dataStart, casacore::Vector<casacore::Int>& dataStep,
112 : casacore::MRadialVelocity& mDataStart,
113 : casacore::MRadialVelocity& mDataStep);
114 : //Utility function to check coordinate match with existing image
115 :
116 : virtual casacore::Bool checkCoord(const casacore::CoordinateSystem& coordsys,
117 : const casacore::String& imageName);
118 :
119 : virtual void setImageParam(casacore::Int& nx, casacore::Int& ny, casacore::Int& npol, casacore::Int& nchan);
120 :
121 : //VisSet and resort
122 : virtual void makeVisSet(casacore::MeasurementSet& ms,
123 : casacore::Bool compress=false, casacore::Bool mosaicOrder=false);
124 : //Just to create the SORTED_TABLE if one can
125 : //virtual void makeVisSet(casacore::MeasurementSet& ms,
126 : // casacore::Bool compress=false, casacore::Bool mosaicOrder=false);
127 :
128 : virtual void writeHistory(casacore::LogIO& os);
129 :
130 : virtual void writeCommand(casacore::LogIO& os);
131 :
132 : //make an empty image
133 : casacore::Bool makeEmptyImage(casacore::CoordinateSystem& imageCoord, casacore::String& name, casacore::Int fieldID=0);
134 :
135 : //Functions to make Primary beams
136 : casacore::Bool makePBImage(casacore::ImageInterface<casacore::Float>& pbImage,
137 : casacore::Bool useSymmetricBeam=true);
138 : casacore::Bool makePBImage(const casacore::CoordinateSystem& imageCoord,
139 : const casacore::String& telescopeName, const casacore::String& diskPBName,
140 : casacore::Bool useSymmetricBeam=true, casacore::Double dishdiam=-1.0);
141 :
142 : casacore::Bool makePBImage(const casacore::CoordinateSystem& imageCoord,
143 : const casacore::Table& vpTable, const casacore::String& diskPBName);
144 :
145 : casacore::Bool makePBImage(const casacore::Table& vpTable, casacore::ImageInterface<casacore::Float>& pbImage);
146 :
147 : casacore::Bool makePBImage(const casacore::CoordinateSystem& imageCoord, PBMath& pbMath, const casacore::String& diskPBName);
148 :
149 : casacore::Bool makePBImage(PBMath& pbMath, casacore::ImageInterface<casacore::Float>& pbImage);
150 :
151 : void setObsInfo(casacore::ObsInfo& obsinfo);
152 : casacore::ObsInfo& latestObsInfo();
153 : // Close the current ms, and replace it with the supplied ms.
154 : // Optionally compress the attached calibration data
155 : // columns if they are created here.
156 : casacore::Bool open(casacore::MeasurementSet &thems, casacore::Bool compress=false, casacore::Bool useModel=false);
157 :
158 : // Flush the ms to disk and detach from the ms file. All function
159 : // calls after this will be a no-op.
160 : casacore::Bool close();
161 :
162 : // Return the name of the MeasurementSet
163 : casacore::String name() const;
164 :
165 : // The following setup methods define the state of the imager.
166 : // <group>
167 : // Set image construction parameters
168 : virtual casacore::Bool setimage(const casacore::Int nx, const casacore::Int ny,
169 : const casacore::Quantity& cellx, const casacore::Quantity& celly,
170 : const casacore::String& stokes,
171 : casacore::Bool doShift,
172 : const casacore::MDirection& phaseCenter,
173 : const casacore::Quantity& shiftx, const casacore::Quantity& shifty,
174 : const casacore::String& mode, const casacore::Int nchan,
175 : const casacore::Int start, const casacore::Int step,
176 : const casacore::MRadialVelocity& mStart, const casacore::MRadialVelocity& mStep,
177 : const casacore::Vector<casacore::Int>& spectralwindowids, const casacore::Int fieldid,
178 : const casacore::Int facets, const casacore::Quantity& distance);
179 :
180 : virtual casacore::Bool defineImage(const casacore::Int nx, const casacore::Int ny,
181 : const casacore::Quantity& cellx, const casacore::Quantity& celly,
182 : const casacore::String& stokes,
183 : const casacore::MDirection& phaseCenter,
184 : const casacore::Int fieldid,
185 : const casacore::String& mode, const casacore::Int nchan,
186 : const casacore::Int start, const casacore::Int step,
187 : const casacore::MFrequency& mFreqStart,
188 : const casacore::MRadialVelocity& mStart,
189 : const casacore::Quantity& qStep,
190 : const casacore::Vector<casacore::Int>& spectralwindowids,
191 : const casacore::Int facets=1,
192 : const casacore::Quantity& restFreq=casacore::Quantity(0,"Hz"),
193 : const casacore::MFrequency::Types& mFreqFrame=casacore::MFrequency::LSRK,
194 : const casacore::Quantity& distance=casacore::Quantity(0,"m"),
195 : const casacore::Bool trackSource=false, const casacore::MDirection&
196 : trackDir=casacore::MDirection(casacore::Quantity(0.0, "deg"),
197 : casacore::Quantity(90.0, "deg")),
198 : const casacore::String& projection=casacore::String("SIN"));
199 : // Set the data selection parameters
200 :
201 : // The parameters useModelcol and readonly is dummy here
202 : //as they are useful for the ImagerMultiMS version only
203 : virtual casacore::Bool setDataPerMS(const casacore::String& msname, const casacore::String& mode,
204 : const casacore::Vector<casacore::Int>& nchan,
205 : const casacore::Vector<casacore::Int>& start,
206 : const casacore::Vector<casacore::Int>& step,
207 : const casacore::Vector<casacore::Int>& spectralwindowids,
208 : const casacore::Vector<casacore::Int>& fieldid,
209 : const casacore::String& msSelect="",
210 : const casacore::String& timerng="",
211 : const casacore::String& fieldnames="",
212 : const casacore::Vector<casacore::Int>& antIndex=casacore::Vector<casacore::Int>(),
213 : const casacore::String& antnames="",
214 : const casacore::String& spwstring="",
215 : const casacore::String& uvdist="",
216 : const casacore::String& scan="",
217 : const casacore::String& intent="",
218 : const casacore::String& obs="",
219 : const casacore::Bool useModelCol=false,
220 : const casacore::Bool readonly=false);
221 :
222 : // Select some data.
223 : // Sets nullSelect_p and returns !nullSelect_p.
224 : // be_calm: lowers the logging level of some messages if true.
225 : casacore::Bool setdata(const casacore::String& mode, const casacore::Vector<casacore::Int>& nchan,
226 : const casacore::Vector<casacore::Int>& start,
227 : const casacore::Vector<casacore::Int>& step, const casacore::MRadialVelocity& mStart,
228 : const casacore::MRadialVelocity& mStep,
229 : const casacore::Vector<casacore::Int>& spectralwindowids,
230 : const casacore::Vector<casacore::Int>& fieldid,
231 : const casacore::String& msSelect="",
232 : const casacore::String& timerng="",
233 : const casacore::String& fieldnames="",
234 : const casacore::Vector<casacore::Int>& antIndex=casacore::Vector<casacore::Int>(),
235 : const casacore::String& antnames="",
236 : const casacore::String& spwstring="",
237 : const casacore::String& uvdist="",
238 : const casacore::String& scan="",
239 : const casacore::String& intent="",
240 : const casacore::String& obs="",
241 : const casacore::Bool usemodelCol=false,
242 : const casacore::Bool be_calm=false);
243 :
244 : // Set the processing options
245 : casacore::Bool setoptions(const casacore::String& ftmachine, const casacore::Long cache, const casacore::Int tile,
246 : const casacore::String& gridfunction, const casacore::MPosition& mLocation,
247 : const casacore::Float padding,
248 : const casacore::Int wprojplanes=-1,
249 : const casacore::String& epJTableName="",
250 : const casacore::Bool applyPointingOffsets=true,
251 : const casacore::Bool doPointingCorrection=true,
252 : const casacore::String& cfCacheDirName="",
253 : const casacore::Float& rotpastep=5.0,
254 : const casacore::Float& computepastep=360.0,
255 : const casacore::Float& pbLimit=5.0e-2,
256 : const casacore::String& freqinterpmethod="linear",
257 : const casacore::Int imageTileSizeInPix=0,
258 : const casacore::Bool singleprecisiononly=false,
259 : const casacore::Int numthreads=-1,
260 : const casacore::Bool psTermOn=true,
261 : const casacore::Bool aTermOn=true,
262 : const casacore::Bool mTermOn=false,
263 : const casacore::Bool wbAWP=false,
264 : const casacore::Bool conjBeams=true);
265 :
266 : // Set the single dish processing options
267 : casacore::Bool setsdoptions(const casacore::Float scale, const casacore::Float weight,
268 : const casacore::Int convsupport=-1, casacore::String pointingColToUse="DIRECTION",
269 : const casacore::Quantity truncate=casacore::Quantity(),
270 : const casacore::Quantity gwidth=casacore::Quantity(),
271 : const casacore::Quantity jwidth=casacore::Quantity(),
272 : const casacore::Float minweight=0.,
273 : const casacore::Bool clipminmax=false,
274 : const casacore::Bool enablecache=false,
275 : const casacore::String & convertfirst="NEVER");
276 :
277 : // Set the voltage pattern
278 : casacore::Bool setvp(const casacore::Bool dovp,
279 : const casacore::Bool defaultVP,
280 : const casacore::String& vpTable,
281 : const casacore::Bool doSquint,
282 : const casacore::Quantity &parAngleInc,
283 : const casacore::Quantity &skyPosThreshold,
284 : casacore::String defaultTel="",
285 : const casacore::Bool verbose=true);
286 :
287 : // Set the scales to be searched in Multi Scale clean
288 : casacore::Bool setscales(const casacore::String& scaleMethod, // "nscales" or "uservector"
289 : const casacore::Int inscales,
290 : const casacore::Vector<casacore::Float>& userScaleSizes);
291 : // set bias
292 : casacore::Bool setSmallScaleBias(const casacore::Float inbias);
293 :
294 : // Set the number of taylor series terms in the expansion of the
295 : // image as a function of frequency.
296 : casacore::Bool settaylorterms(const casacore::Int intaylor,
297 : const casacore::Double inreffreq);
298 :
299 : // </group>
300 :
301 : // Advise on suitable values
302 : casacore::Bool advise(const casacore::Bool takeAdvice, const casacore::Float amplitudeloss,
303 : const casacore::Quantity& fieldOfView,
304 : casacore::Quantity& cell, casacore::Int& npixels, casacore::Int& facets,
305 : casacore::MDirection& phaseCenter);
306 :
307 : // Output a summary of the state of the object
308 : casacore::Bool summary();
309 :
310 : // Return the state of the object as a string
311 : casacore::String state();
312 :
313 : // Return the # of visibilities accessible to *rvi, optionally excluding
314 : // flagged ones (if unflagged_only is true) and/or ones without imaging
315 : // weights (if must_have_imwt is true).
316 : casacore::uInt count_visibilities(ROVisibilityIterator *rvi,
317 : const casacore::Bool unflagged_only, const casacore::Bool must_have_imwt);
318 :
319 : // Return the image coordinates
320 : casacore::Bool imagecoordinates(casacore::CoordinateSystem& coordInfo, const casacore::Bool verbose=true);
321 : // new version
322 : casacore::Bool imagecoordinates2(casacore::CoordinateSystem& coordInfo, const casacore::Bool verbose=true);
323 :
324 : // Return the image shape
325 : casacore::IPosition imageshape() const;
326 :
327 : // Weight the MeasurementSet
328 : //For some time of weighting briggs/uniform ...one can do it on a per field basis to calculate
329 : //weight density distribution. If that is what is wanted multiField should be set to true
330 : //multifield is inoperative for natural, radial weighting
331 : casacore::Bool weight(const casacore::String& algorithm, const casacore::String& rmode,
332 : const casacore::Quantity& noise, const casacore::Double robust,
333 : const casacore::Quantity& fieldofview, const casacore::Int npixels, const casacore::Bool multiField=false);
334 :
335 : // Filter the MeasurementSet
336 : casacore::Bool filter(const casacore::String& type, const casacore::Quantity& bmaj, const casacore::Quantity& bmin,
337 : const casacore::Quantity& bpa);
338 :
339 : // Apply a uvrange
340 : casacore::Bool uvrange(const casacore::Double& uvmin, const casacore::Double& uvmax);
341 :
342 : // Sensitivity
343 : casacore::Bool sensitivity(casacore::Quantity& pointsourcesens, casacore::Double& relativesens, casacore::Double& sumwt,
344 : casacore::Double& effectiveBandwidth, casacore::Double& effectiveIntegration, casacore::Int& mBaselines,
345 : casacore::Matrix<casacore::Int>& mssChanSel,
346 : casacore::Vector<casacore::Vector<casacore::Int> >& nData,
347 : casacore::Vector<casacore::Vector<casacore::Double> >& sumwtChan,
348 : casacore::Vector<casacore::Vector<casacore::Double> >& sumwtsqChan,
349 : casacore::Vector<casacore::Vector<casacore::Double> >& sumInverseVarianceChan);
350 :
351 : // Apparent sensitivity calculator
352 : // Accounts for gridding weights and assumes casacore::MS weights have correct units
353 : casacore::Bool apparentSensitivity(casacore::Double& effSensitivity,
354 : casacore::Double& relToNat);
355 :
356 : // Make plain image + keep the complex image as well if complexImageName != "".
357 : casacore::Bool makeimage(const casacore::String& type, const casacore::String& imageName,
358 : const casacore::String& complexImageName="", const casacore::Bool verbose=true);
359 :
360 : // Fill in a region of a mask
361 : casacore::Bool boxmask(const casacore::String& mask, const casacore::Vector<casacore::Int>& blc,
362 : const casacore::Vector<casacore::Int>& trc,const casacore::Float value);
363 :
364 : //Make a region either from record or array of blc trc
365 : //(casacore::Matrix(nboxes,4)) into a mask image
366 : //value is the value of the mask pixels
367 : //circular masks has form casacore::Matrix(ncircles,3)
368 : //where the 3 values on a row are radius, x, y pixel values
369 : casacore::Bool regionmask(const casacore::String& maskimage, casacore::Record* imageRegRec,
370 : casacore::Matrix<casacore::Quantity>& blctrcs, casacore::Matrix<casacore::Float>& circles,
371 : const casacore::Float& value=1.0);
372 :
373 : static casacore::Bool regionToImageMask(const casacore::String& maskimage, casacore::Record* imageRegRec,
374 : casacore::Matrix<casacore::Quantity>& blctrcs,
375 : casacore::Matrix<casacore::Float>& circles,
376 : const casacore::Float& value=1.0);
377 : // Clip on casacore::Stokes I
378 : casacore::Bool clipimage(const casacore::String& image, const casacore::Quantity& threshold);
379 :
380 : // Make a mask image
381 : static casacore::Bool mask(const casacore::String& mask, const casacore::String& imageName,
382 : const casacore::Quantity& threshold);
383 :
384 : // Restore
385 : casacore::Bool restore(const casacore::Vector<casacore::String>& model, const casacore::String& complist,
386 : const casacore::Vector<casacore::String>& image, const casacore::Vector<casacore::String>& residual);
387 :
388 : // similar to restore except this is to be called if you fiddle with the model and complist
389 : // outside of this object (say you clip stuff etc) ...keep the sm_p and se_p state but just calculate new residuals and
390 : // restored images. Will throw an exception is se_p or sm_p is not valid (i.e you should have used clean, mem etc before hand).
391 : casacore::Bool updateresidual(const casacore::Vector<casacore::String>& model, const casacore::String& complist,
392 : const casacore::Vector<casacore::String>& image, const casacore::Vector<casacore::String>& residual);
393 :
394 : // Setbeam
395 : casacore::Bool setbeam(const casacore::ImageBeamSet& beam);
396 :
397 : // Residual
398 : casacore::Bool residual(const casacore::Vector<casacore::String>& model, const casacore::String& complist,
399 : const casacore::Vector<casacore::String>& image);
400 :
401 : // Approximate PSF
402 : casacore::Bool approximatepsf(const casacore::String& psf);
403 :
404 : // Smooth
405 : casacore::Bool smooth(const casacore::Vector<casacore::String>& model,
406 : const casacore::Vector<casacore::String>& image, casacore::Bool usefit,
407 : casacore::ImageBeamSet& mbeam,
408 : casacore::Bool normalizeVolume);
409 :
410 : // Clean algorithm
411 : casacore::Record clean(const casacore::String& algorithm,
412 : const casacore::Int niter,
413 : const casacore::Float gain,
414 : const casacore::Quantity& threshold,
415 : const casacore::Bool displayProgress,
416 : const casacore::Vector<casacore::String>& model, const casacore::Vector<casacore::Bool>& fixed,
417 : const casacore::String& complist,
418 : const casacore::Vector<casacore::String>& mask,
419 : const casacore::Vector<casacore::String>& restored,
420 : const casacore::Vector<casacore::String>& residual,
421 : const casacore::Vector<casacore::String>& psf=casacore::Vector<casacore::String>(0),
422 : const casacore::Bool firstrun=true);
423 :
424 : casacore::Record iClean(const casacore::String& algorithm,
425 : const casacore::Int niter,
426 : const casacore::Double gain,
427 : //const casacore::String& threshold,
428 : const casacore::Quantity& threshold,
429 : const casacore::Bool displayprogress,
430 : const casacore::Vector<casacore::String>& model,
431 : const casacore::Vector<casacore::Bool>& keepfixed, const casacore::String& complist,
432 : const casacore::Vector<casacore::String>& mask,
433 : const casacore::Vector<casacore::String>& image,
434 : const casacore::Vector<casacore::String>& residual,
435 : const casacore::Vector<casacore::String>& psfnames,
436 : const casacore::Bool interactive, const casacore::Int npercycle,
437 : const casacore::String& masktemplate);
438 :
439 : // MEM algorithm
440 : casacore::Bool mem(const casacore::String& algorithm,
441 : const casacore::Int niter, const casacore::Quantity& sigma,
442 : const casacore::Quantity& targetflux,
443 : const casacore::Bool constrainflux,
444 : const casacore::Bool displayProgress,
445 : const casacore::Vector<casacore::String>& model, const casacore::Vector<casacore::Bool>& fixed,
446 : const casacore::String& complist,
447 : const casacore::Vector<casacore::String>& prior,
448 : const casacore::Vector<casacore::String>& mask,
449 : const casacore::Vector<casacore::String>& restored,
450 : const casacore::Vector<casacore::String>& residual);
451 :
452 : // NNLS algorithm
453 : casacore::Bool nnls(const casacore::String& algorithm, const casacore::Int niter, const casacore::Float tolerance,
454 : const casacore::Vector<casacore::String>& model, const casacore::Vector<casacore::Bool>& fixed,
455 : const casacore::String& complist,
456 : const casacore::Vector<casacore::String>& fluxMask, const casacore::Vector<casacore::String>& dataMask,
457 : const casacore::Vector<casacore::String>& restored,
458 : const casacore::Vector<casacore::String>& residual);
459 :
460 : // Multi-field control parameters
461 : //flat noise is the parameter that control the search of clean components
462 : //in a flat noise image or an optimum beam^2 image
463 : casacore::Bool setmfcontrol(const casacore::Float cyclefactor,
464 : const casacore::Float cyclespeedup,
465 : const casacore::Float cyclemaxpsffraction,
466 : const casacore::Int stoplargenegatives,
467 : const casacore::Int stoppointmode,
468 : const casacore::String& scaleType,
469 : const casacore::Float minPB,
470 : const casacore::Float constPB,
471 : const casacore::Vector<casacore::String>& fluxscale,
472 : const casacore::Bool flatnoise=true);
473 :
474 : // Feathering algorithm
475 : casacore::Bool feather(const casacore::String& image,
476 : const casacore::String& highres,
477 : const casacore::String& lowres,
478 : const casacore::String& lowpsf, const casacore::Float dishDiam=-1.0, const casacore::Bool lowPassFilterSD=false);
479 :
480 : // Apply or correct for Primary Beam or Voltage Pattern
481 : casacore::Bool pb(const casacore::String& inimage,
482 : const casacore::String& outimage,
483 : const casacore::String& incomps,
484 : const casacore::String& outcomps,
485 : const casacore::String& operation,
486 : const casacore::MDirection& pointngCenter,
487 : const casacore::Quantity& pa,
488 : const casacore::String& pborvp);
489 :
490 : // Make a linear mosaic of several images
491 : casacore::Bool linearmosaic(const casacore::String& mosaic,
492 : const casacore::String& fluxscale,
493 : const casacore::String& sensitivity,
494 : const casacore::Vector<casacore::String>& images,
495 : const casacore::Vector<casacore::Int>& fieldids);
496 :
497 : // Fourier transform the model and componentlist.
498 : // Returns its nominal success value.
499 : ///For moving time variable phasecenters in field table
500 : //phasecentertime is usually the time in the data unless
501 : //you want to use a specific time which you can set here
502 : casacore::Bool ft(const casacore::Vector<casacore::String>& model, const casacore::String& complist,
503 : casacore::Bool incremental=false, const casacore::Double phasecenterTime=-1.0);
504 :
505 : // Compute the model visibility using specified source flux densities
506 : casacore::Bool setjy(const casacore::Int fieldid, const casacore::Int spectralwindowid,
507 : const casacore::Vector<casacore::Double>& fluxDensity, const casacore::String& standard);
508 : casacore::Bool setjy(const casacore::Vector<casacore::Int>& fieldid, const casacore::Vector<casacore::Int>& spectralwindowid,
509 : const casacore::String& fieldnames, const casacore::String& spwstring,
510 : const casacore::Vector<casacore::Double>& fluxDensity, const casacore::String& standard);
511 :
512 : //Setjy with model image. If chanDep=true then the scaling is calculated on a
513 : //per channel basis for the model image...otherwise the whole spw gets the
514 : //same flux density.
515 : //casacore::Bool setjy(const casacore::Vector<casacore::Int>& fieldid,
516 : casacore::Record setjy(const casacore::Vector<casacore::Int>& fieldid,
517 : const casacore::Vector<casacore::Int>& spectralwindowid,
518 : const casacore::String& fieldnames, const casacore::String& spwstring,
519 : const casacore::String& model,
520 : const casacore::Vector<casacore::Double>& fluxDensity, const casacore::String& standard,
521 : const casacore::Bool chanDep=false, //const casacore::Double spix=0.0,
522 : const casacore::Vector<casacore::Double>& spix=casacore::Vector<casacore::Double>(),
523 : const casacore::MFrequency& reffreq=casacore::MFrequency(casacore::Quantity(1.0, "GHz"),
524 : casacore::MFrequency::LSRK),
525 : const casacore::Vector<casacore::Double>& pipars=casacore::Vector<casacore::Double>(),
526 : const casacore::Vector<casacore::Double>& papars=casacore::Vector<casacore::Double>(),
527 : const casacore::Double& rotMeas=0.0,
528 : const casacore::String& timerange="", const casacore::String& scanstr="",
529 : const casacore::String& intentstr="", const casacore::String& obsidstr="",
530 : const casacore::String& interpolation="nearest");
531 :
532 : // Make an empty image
533 : casacore::Bool make(const casacore::String& model);
534 :
535 : // make a model from a SD image.
536 : // This model then can be used as initial clean model to include the
537 : // shorter spacing.
538 : casacore::Bool makemodelfromsd(const casacore::String& sdImage, const casacore::String& modelimage,
539 : const casacore::String& lowPSF,
540 : casacore::String& maskImage);
541 :
542 : // Write a component list to disk, starting with prefix, using a setjy
543 : // standard, and return the name of the list.
544 : casacore::String make_comp(const casacore::String& objName, const casacore::String& standard,
545 : const casacore::MEpoch& mtime, const casacore::Vector<casacore::MFrequency>& freqv,
546 : const casacore::String& prefix);
547 :
548 : // Clone an image
549 : static casacore::Bool clone(const casacore::String& imageName, const casacore::String& newImageName);
550 :
551 : // Fit the psf
552 : casacore::Bool fitpsf(const casacore::String& psf, casacore::ImageBeamSet& mbeam);
553 :
554 : // Correct the visibility data (OBSERVED->CORRECTED)
555 : casacore::Bool correct(const casacore::Bool doparallactic, const casacore::Quantity& t);
556 :
557 : // Plot the uv plane
558 : casacore::Bool plotuv(const casacore::Bool rotate);
559 :
560 : // Plot the visibilities
561 : casacore::Bool plotvis(const casacore::String& type, const casacore::Int increment);
562 :
563 : // Plot the weights
564 : casacore::Bool plotweights(const casacore::Bool gridded, const casacore::Int increment);
565 :
566 : // Plot a summary
567 : casacore::Bool plotsummary();
568 :
569 : // Clip visibilities
570 : casacore::Bool clipvis(const casacore::Quantity& threshold);
571 :
572 :
573 : //Check if can proceed with this object
574 : casacore::Bool valid() const;
575 :
576 :
577 : //Interactive mask drawing
578 : //forceReload..forces the viewer to dump previous image that is being displayed
579 : casacore::Int interactivemask(const casacore::String& imagename, const casacore::String& maskname,
580 : casacore::Int& niter, casacore::Int& ncycles, casacore::String& threshold, const casacore::Bool forceReload=false);
581 :
582 :
583 : //helper function to copy a mask from one image to another
584 :
585 : static casacore::Bool copyMask(casacore::ImageInterface<casacore::Float>& out, const casacore::ImageInterface<casacore::Float>& in, casacore::String maskname="mask0", casacore::Bool setdefault=true);
586 :
587 :
588 : // Supports the "[] or -1 => everything" convention using the rule:
589 : // If v is empty or only has 1 element, and it is < 0,
590 : // replace v with 0, 1, ..., nelem - 1.
591 : // Returns whether or not it modified v.
592 : // If so, v is modified in place.
593 : static casacore::Bool expand_blank_sel(casacore::Vector<casacore::Int>& v, const casacore::uInt nelem);
594 :
595 : //spectral gridding calculation for output images (use SubMS::calcChanFreqs)
596 : casacore::Bool calcImFreqs(casacore::Vector<casacore::Double>& imfreqs, casacore::Vector<casacore::Double>& imfreqres,
597 : const casacore::MFrequency::Types& oldRefFrame,
598 : const casacore::MEpoch& obsEpoch, const casacore::MPosition& obsPosition,
599 : const casacore::Double& restFreq);
600 :
601 : // Advise the chanselection needed for the frequency range or
602 : // give the frequency range for a give spwselection if getFreqRange==true
603 : // if the parameter msname is used then the MSs associated associated with
604 : // this object (that have been either 'open'ed or 'selectvis'ed) are ignored
605 : // In this mode it is a helper function to the general world ...no need to
606 : // open or selectvis. You need to specify the field_id for which this calculation is
607 : // being done for in the helper mode.
608 : // If you have already set casacore::MS's and selected data and msname="" then
609 : // the calulation is done for the field(s) selected in selectvis.
610 : // getFreqRange=true then the freqrange in the frame and spwselection you choose is
611 : // returned in freqStart and freqEnd (in the case of msname="" then it is for the fields
612 : //and spw you have chosen in selectvis).
613 : casacore::Bool adviseChanSelex(casacore::Double& freqStart, casacore::Double& freqEnd,
614 : const casacore::Double& freqStep, const casacore::MFrequency::Types& freqframe,
615 : casacore::Vector< casacore::Vector<casacore::Int> >& spw, casacore::Vector< casacore::Vector<casacore::Int> >& start,
616 : casacore::Vector< casacore::Vector<casacore::Int> >& nchan, const casacore::String& msname="",
617 : const casacore::Int fieldid=0, const casacore::Bool getFreqRange=false,
618 : const casacore::String spwselection="");
619 :
620 :
621 : //These are utility functions when weights from different imager instances
622 : //need to reconciled in parallel gridding by different instances of imagers
623 : //for example.
624 : // when type is "imaging"
625 : // getweightGrid will get the weight density for uniform style imaging weight
626 : // the casacore::Block elements are for different fields if independent field weighting
627 : // was done.
628 : // when type is "ftweight"..then a casacore::Vector of string is expected in weightimage
629 : // which is of the same length as the number of models put in clean etc
630 :
631 : casacore::Bool getWeightGrid(casacore::Block<casacore::Matrix<casacore::Float> >&weightgrid, const casacore::String& type, const casacore::Vector<casacore::String>& weightImagenames=casacore::Vector<casacore::String>());
632 : casacore::Bool setWeightGrid(const casacore::Block<casacore::Matrix<casacore::Float> >& weightgrid, const casacore::String& type);
633 : casacore::String dQuantitytoString(const casacore::Quantity& dq);
634 :
635 : // Automatic evaluation of map extent for given visibility.
636 : //
637 : //
638 : // @param[in] referenceFrame reference direction frame
639 : // @param[in] movingSource name of moving source
640 : // @param[in] pointingColumn pointing column to use
641 : // @param[out] center center of the map
642 : // @param[out] blc bottom left corner of the map
643 : // @param[out] trc top right corner of the map
644 : // @param[out] extent map extent
645 : //
646 : // @return
647 : virtual casacore::Bool mapExtent(const casacore::String &referenceFrame, const casacore::String &movingSource,
648 : const casacore::String &pointingColumn, casacore::Vector<casacore::Double> ¢er, casacore::Vector<casacore::Double> &blc,
649 : casacore::Vector<casacore::Double> &trc, casacore::Vector<casacore::Double> &extent);
650 :
651 : // Caluculate sampling interval of Raster map
652 : //
653 : //
654 : // @param[in] referenceFrame reference direction frame
655 : // @param[in] movingSource name of moving source
656 : // @param[in] pointingColumn pointing column to use
657 : // @param[in] antenna antenna selection string, e.g., '0&&&'
658 : // @param[out] sampling the sampling interval in scan direction and orthogonal direction
659 : // @param[out] positionAngle position angle of scan
660 : //
661 : // @return
662 : virtual casacore::Bool pointingSampling(const casacore::String &referenceFrame, const casacore::String &movingSource,
663 : const casacore::String &pointingColumn, const casacore::String &antenna,
664 : casacore::Quantum<casacore::Vector<casacore::Double>> &sampling, casacore::Quantity &positionAngle);
665 :
666 : //Helper function to transfer history table to a logger holder
667 : //which can be stored in images
668 : static void transferHistory(casacore::LoggerHolder& imageLog, casacore::MSHistoryColumns& msHis);
669 :
670 : protected:
671 :
672 : casacore::CountedPtr<casacore::MeasurementSet> ms_p;
673 : casacore::CountedPtr<casacore::MSHistoryHandler> hist_p;
674 : casacore::Table antab_p;
675 : casacore::Table datadesctab_p;
676 : casacore::Table feedtab_p;
677 : casacore::Table fieldtab_p;
678 : casacore::Table obstab_p;
679 : casacore::Table pointingtab_p;
680 : casacore::Table poltab_p;
681 : casacore::Table proctab_p;
682 : casacore::Table spwtab_p;
683 : casacore::Table statetab_p;
684 : casacore::Table dopplertab_p;
685 : casacore::Table flagcmdtab_p;
686 : casacore::Table freqoffsettab_p;
687 : casacore::Table historytab_p;
688 : casacore::Table sourcetab_p;
689 : casacore::Table syscaltab_p;
690 : casacore::Table weathertab_p;
691 : casacore::Int lockCounter_p;
692 : casacore::Int nx_p, ny_p, npol_p, nchan_p;
693 : casacore::ObsInfo latestObsInfo_p;
694 : //What should be the tile volume on disk
695 : casacore::Int imageTileVol_p;
696 :
697 :
698 :
699 : casacore::String msname_p;
700 : casacore::CountedPtr<casacore::MeasurementSet> mssel_p;
701 : VisSet *vs_p;
702 : ROVisibilityIterator* rvi_p;
703 : VisibilityIterator* wvi_p;
704 : FTMachine *ft_p;
705 : ComponentFTMachine *cft_p;
706 : SkyEquation* se_p;
707 : CleanImageSkyModel* sm_p;
708 : VPSkyJones* vp_p;
709 : VPSkyJones* gvp_p;
710 :
711 : casacore::Bool setimaged_p, nullSelect_p;
712 : casacore::Bool redoSkyModel_p; // if clean is run multiply ..use this to check
713 : // if setimage was changed hence redo the skyModel.
714 : casacore::Float rotPAStep_p, computePAStep_p, pbLimit_p;
715 : casacore::Int facets_p;
716 : casacore::Int wprojPlanes_p;
717 : casacore::Quantity mcellx_p, mcelly_p;
718 : casacore::String stokes_p;
719 : casacore::String dataMode_p;
720 : casacore::String imageMode_p; // channel, (optical)velocity, mfs, or frequency
721 : casacore::Vector<casacore::Int> dataNchan_p;
722 : casacore::Int imageNchan_p;
723 : casacore::Vector<casacore::Int> dataStart_p, dataStep_p;
724 : casacore::Int imageStart_p, imageStep_p;
725 : casacore::MRadialVelocity mDataStart_p, mImageStart_p;
726 : casacore::MRadialVelocity mDataStep_p, mImageStep_p;
727 : casacore::MFrequency mfImageStart_p, mfImageStep_p;
728 : casacore::MFrequency::Types freqFrame_p;
729 : casacore::MDirection phaseCenter_p;
730 : casacore::Quantity restFreq_p;
731 : casacore::Quantity distance_p;
732 : casacore::Bool doShift_p;
733 : casacore::Quantity shiftx_p;
734 : casacore::Quantity shifty_p;
735 : casacore::String ftmachine_p, gridfunction_p;
736 : casacore::Bool wfGridding_p;
737 : casacore::Long cache_p;
738 : casacore::Int tile_p;
739 : casacore::MPosition mLocation_p;
740 : casacore::Bool doVP_p;
741 : casacore::ImageBeamSet beam_p;
742 : casacore::Bool beamValid_p;
743 : casacore::Float padding_p;
744 : casacore::Float sdScale_p;
745 : casacore::Float sdWeight_p;
746 : casacore::Int sdConvSupport_p;
747 :
748 : casacore::Quantity qtruncate_p;
749 : casacore::Quantity qgwidth_p;
750 : casacore::Quantity qjwidth_p;
751 :
752 : casacore::Float minWeight_p;
753 :
754 : casacore::Bool clipminmax_p;
755 :
756 : casacore::Bool enablecache_p;
757 :
758 : casacore::String convertfirst_p;
759 :
760 : // special mf control parms, etc
761 : casacore::Float cyclefactor_p;
762 : casacore::Float cyclespeedup_p;
763 : casacore::Float cyclemaxpsffraction_p;
764 : casacore::Int stoplargenegatives_p;
765 : casacore::Int stoppointmode_p;
766 : casacore::Vector<casacore::String> fluxscale_p;
767 : casacore::String scaleType_p; // type of image-plane scaling: NONE, SAULT
768 : casacore::Float minPB_p; // minimum value of generalized-PB pattern
769 : casacore::Float constPB_p; // above this level, constant flux-scale
770 :
771 : casacore::Vector<casacore::Int> spectralwindowids_p;
772 : casacore::Int fieldid_p;
773 :
774 : casacore::Vector<casacore::Int> dataspectralwindowids_p;
775 : casacore::Vector<casacore::Int> datadescids_p;
776 : casacore::Vector<casacore::Int> datafieldids_p;
777 : //TT
778 : casacore::Cube<casacore::Int> spwchansels_p;
779 : casacore::Matrix<casacore::Double> freqrange_p;
780 : casacore::Matrix<casacore::Double> mssFreqSel_p;
781 : casacore::Matrix<casacore::Int> mssChanSel_p;
782 :
783 : casacore::Int numMS_p;
784 :
785 : casacore::String telescope_p;
786 : casacore::String vpTableStr_p; // description of voltage patterns for various telescopes
787 : // in the MS
788 : casacore::Quantity parAngleInc_p;
789 : casacore::Quantity skyPosThreshold_p;
790 : BeamSquint::SquintType squintType_p;
791 : casacore::Bool doDefaultVP_p; // make default VPs, rather than reading in a vpTable
792 :
793 :
794 : casacore::Bool doMultiFields_p; // Do multiple fields?
795 : casacore::Bool multiFields_p; // multiple fields have been specified in setdata
796 :
797 : casacore::Bool doWideBand_p; // Do Multi Frequency Synthesis Imaging
798 : casacore::String freqInterpMethod_p; //frequency interpolation mode
799 :
800 : casacore::Bool flatnoise_p;
801 :
802 : // Set the defaults
803 : void defaults();
804 :
805 : // check if it is dettahced from ms.
806 : casacore::Bool detached() const;
807 :
808 : // Create the FTMachines when necessary or when the control parameters
809 : // have changed.
810 : virtual casacore::Bool createFTMachine();
811 :
812 : void openSubTable (const casacore::Table & otherTable, casacore::Table & table, const casacore::TableLock & tableLock);
813 :
814 : casacore::Bool removeTable(const casacore::String& tablename);
815 : casacore::Bool updateSkyModel(const casacore::Vector<casacore::String>& model,
816 : const casacore::String complist);
817 : casacore::Bool createSkyEquation(const casacore::String complist="");
818 : casacore::Bool createSkyEquation(const casacore::Vector<casacore::String>& image,
819 : const casacore::Vector<casacore::Bool>& fixed,
820 : const casacore::String complist="");
821 : casacore::Bool createSkyEquation(const casacore::Vector<casacore::String>& image,
822 : const casacore::String complist="");
823 : casacore::Bool createSkyEquation(const casacore::Vector<casacore::String>& image,
824 : const casacore::Vector<casacore::Bool>& fixed,
825 : const casacore::Vector<casacore::String>& mask,
826 : const casacore::String complist="");
827 : casacore::Bool createSkyEquation(const casacore::Vector<casacore::String>& image,
828 : const casacore::Vector<casacore::Bool>& fixed,
829 : const casacore::Vector<casacore::String>& mask,
830 : const casacore::Vector<casacore::String>& fluxMask,
831 : const casacore::String complist="");
832 : ATerm* createTelescopeATerm(casacore::MeasurementSet& ms, const casacore::Bool& isATermOn=true);
833 : void destroySkyEquation();
834 :
835 : //add residual to the private vars or create residual images
836 : casacore::Bool addResiduals(const casacore::Vector<casacore::String>& residual);
837 : // Add the residuals to the SkyEquation
838 : casacore::Bool addResidualsToSkyEquation(const casacore::Vector<casacore::String>& residual);
839 :
840 : // Add or replace the masks
841 : casacore::Bool addMasksToSkyEquation(const casacore::Vector<casacore::String>& mask, const casacore::Vector<casacore::Bool>& fixed=casacore::Vector<casacore::Bool>(0));
842 :
843 : // Get the rest frequency ..returns 1 element in restfreq
844 : // if user specified or try to get the info from the SOURCE table
845 : casacore::Bool getRestFreq(casacore::Vector<casacore::Double>& restFreq, const casacore::Int& spw);
846 :
847 : casacore::Bool restoreImages(const casacore::Vector<casacore::String>& restored, casacore::Bool modresiduals=true);
848 :
849 : // names of flux scale images
850 : casacore::Bool writeFluxScales(const casacore::Vector<casacore::String>& fluxScaleNames);
851 :
852 : // Helper functions to hide some setjy code.
853 : casacore::Unit sjy_setup_arrs(casacore::Vector<casacore::Vector<Flux<casacore::Double> > >& returnFluxes,
854 : casacore::Vector<casacore::Vector<Flux<casacore::Double> > >& returnFluxErrs,
855 : casacore::Vector<casacore::Vector<casacore::Double> >& fluxUsed, // mainly for logging purpose
856 : casacore::Vector<casacore::String>& tempCLs,
857 : casacore::Vector<casacore::Vector<casacore::MFrequency> >& mfreqs,
858 : const casacore::MSSpWindowColumns& spwcols, const casacore::uInt nspws,
859 : const casacore::Vector<casacore::Int>& selToRawSpwIds, const casacore::Bool chanDep);
860 : // Returns whether it might have made any visibilities.
861 : casacore::Bool sjy_make_visibilities(casacore::TempImage<casacore::Float> *tmodimage, casacore::LogIO& os,
862 : //casacore::Bool sjy_make_visibilities(casacore::Block<casacore::CountedPtr<casacore::TempImage<casacore::Float> > >& tmodimages, casacore::LogIO& os,
863 : //const casacore::Int rawspwid, const casacore::Int fldid,
864 : // for new one
865 : // const casacore::Vector<casacore::Int>& rawspwids, const casacore::Int fldid,
866 : const casacore::Int rawspwid, const casacore::Int fldid,
867 : const casacore::String& clname, const casacore::String& timerange="",
868 : const casacore::String& scanstr="",
869 : const casacore::String& obsidstr="",
870 : const casacore::String& intentstr="",
871 : const casacore::Vector<casacore::Double>& freqofscale=casacore::Vector<casacore::Double>(0),
872 : const casacore::Vector<casacore::Double>& scale=casacore::Vector<casacore::Double>(0) );
873 : // Concatenate multiple CLs
874 : casacore::Bool sjy_concatComponentLists(casacore::LogIO& os, const casacore::Vector<casacore::String>& tempCLs, const casacore::String& outTempCL);
875 : // Returns whether it found a source.
876 : casacore::Bool sjy_computeFlux(casacore::LogIO& os, FluxStandard& fluxStd,
877 : casacore::Vector<casacore::Vector<Flux<casacore::Double> > >& returnFluxes,
878 : casacore::Vector<casacore::Vector<Flux<casacore::Double> > >& returnFluxErrs,
879 : casacore::Vector<casacore::String>& tempCLs,
880 : //casacore::Vector<casacore::Double>& fluxUsed,
881 : casacore::Vector<casacore::Vector<casacore::Double> >& fluxUsed,
882 : casacore::String& fluxScaleName, casacore::MEpoch& aveEpoch,
883 : const casacore::Vector<casacore::Vector<casacore::MFrequency> >& mfreqs,
884 : const casacore::String& model, const casacore::String& fieldName,
885 : const casacore::MSColumns& msc, const casacore::Int fldid,
886 : const casacore::MDirection& fieldDir, const casacore::Vector<casacore::Int>& selToRawSpwIds,
887 : const casacore::String& standard);
888 :
889 : void sjy_makeComponentList(casacore::LogIO& os, casacore::Vector<casacore::String>& tempCLs,
890 : casacore::Vector<casacore::Vector<Flux<casacore::Double> > >& returnFluxes,
891 : const casacore::Vector<casacore::Double>& fluxUsed,
892 : const casacore::Vector<casacore::Int>& selToRawSpwIds,
893 : const casacore::Vector<casacore::Vector<casacore::MFrequency> >& mfreqs,
894 : const casacore::String& fieldName,
895 : const casacore::MDirection& fieldDir,
896 : const casacore::Vector<casacore::Double>& spix,
897 : const casacore::Vector<casacore::Double>& pipars,
898 : const casacore::Vector<casacore::Double>& papars,
899 : const casacore::Double& rotMeas,
900 : //const casacore::Vector<casacore::Double>& cppars,
901 : const casacore::MFrequency& reffreq,
902 : const casacore::MEpoch& aveEpoch,
903 : const casacore::Int fldId);
904 : //
905 : // Returns NULL if no image is prepared.
906 : casacore::TempImage<casacore::Float>* sjy_prepImage(casacore::LogIO& os, FluxStandard& fluxStd,
907 : casacore::Vector<casacore::Double>& fluxUsed,
908 : casacore::Vector<casacore::Double>& freq,
909 : casacore::Vector<casacore::Double>& scale, const casacore::String& model,
910 : const casacore::MSSpWindowColumns& spwcols,
911 : //const casacore::Int rawspwid, const casacore::Bool chanDep,
912 : const casacore::Vector<casacore::Int> rawspwids, const casacore::Bool chanDep,
913 : const casacore::Vector<casacore::Vector<casacore::MFrequency> >& mfreqs,
914 : //const casacore::uInt selspw, const casacore::String& fieldName,
915 : const casacore::String& fieldName,
916 : const casacore::MDirection& fieldDir, const casacore::Unit& freqUnit,
917 : const casacore::Vector<casacore::Double>& fluxdens,
918 : const casacore::Bool precompute,
919 : //const casacore::Double spix,
920 : const casacore::Vector<casacore::Double>& spix,
921 : const casacore::MFrequency& reffreq,
922 : const casacore::MEpoch& aveEpoch, const casacore::Int fieldId);
923 : // Returns true or throws up.
924 : casacore::Bool sjy_regridCubeChans(casacore::TempImage<casacore::Float>* tmodimage,
925 : casacore::PagedImage<casacore::Float>& modimage, casacore::Int freqAxis);
926 :
927 : // set a radius limit if the model image is one of the known source
928 : // otherwise simply copy modimage to tmodimage
929 : casacore::Bool sjy_setRadiusLimit(casacore::TempImage<casacore::Float>* tmodimage,
930 : casacore::PagedImage<casacore::Float>& modimage, const casacore::String& model,
931 : casacore::DirectionCoordinate& dircsys);
932 :
933 : casacore::Bool sjy_calciflux(const casacore::Vector<casacore::MFrequency>& freqs, const casacore::MFrequency& reffreq,
934 : const casacore::Double refflux, const casacore::Vector<casacore::Double>& vspix, casacore::Vector<casacore::Double>& iflux);
935 :
936 : casacore::Bool sjy_calcquflux(const casacore::Vector<casacore::Double>& pipars, const casacore::Vector<casacore::Double>& papars,
937 : const casacore::Vector<casacore::Double>& iflux, const casacore::Double rotMeas,
938 : const casacore::Vector<casacore::MFrequency>& freqs,
939 : const casacore::MFrequency& reffreq, casacore::Vector<casacore::Double>& qflux,
940 : casacore::Vector<casacore::Double>& uflux);
941 :
942 : casacore::String imageName();
943 :
944 : casacore::Bool pbguts(casacore::ImageInterface<casacore::Float>& in,
945 : casacore::ImageInterface<casacore::Float>& out,
946 : const casacore::MDirection&,
947 : const casacore::Quantity&);
948 :
949 : // Helper func for printing clean's restoring beam to the logger. May find
950 : // the restoring beam as a side effect, so sm_p can't be const.
951 : void printbeam(CleanImageSkyModel *sm_p, casacore::LogIO &os, const casacore::Bool firstrun=true);
952 :
953 : // Helper func for createFTMachine(). Returns phaseCenter_p as a casacore::String,
954 : // *assuming* it is set. It does not check!
955 : casacore::String tangentPoint();
956 :
957 :
958 : casacore::Bool assertDefinedImageParameters() const;
959 : // Virtual methods to set the ImageSkyModel and SkyEquation.
960 : // This allows derived class pimager to set parallelized
961 : // specializations.
962 : //
963 0 : virtual void setWFCleanImageSkyModel()
964 0 : {sm_p = new WFCleanImageSkyModel(facets_p, wfGridding_p); return;};
965 :
966 0 : virtual void setClarkCleanImageSkyModel()
967 0 : {sm_p = new ClarkCleanImageSkyModel(); return;};
968 : virtual void setSkyEquation();
969 :
970 : virtual void savePSF(const casacore::Vector<casacore::String>& psf);
971 :
972 : casacore::String frmtTime(const casacore::Double time);
973 :
974 : //copy imageregion to pixels on image as value given
975 : static casacore::Bool regionToMask(casacore::ImageInterface<casacore::Float>& maskImage, casacore::ImageRegion& imagreg, const casacore::Float& value=1.0);
976 :
977 : //set the mosaic ft machine and right convolution function
978 : virtual void setMosaicFTMachine(casacore::Bool useDoublePrec=false);
979 :
980 : // Makes a component list on disk containing cmp (with fluxval and cspectrum)
981 : // named msname_p.fieldName.spw<spwid>.tempcl and returns the name.
982 : casacore::String makeComponentList(const casacore::String& fieldName, const casacore::Int spwid,
983 : const Flux<casacore::Double>& fluxval,
984 : const ComponentShape& cmp,
985 : const ConstantSpectrum& cspectrum) const;
986 :
987 : casacore::Vector<casacore::Int> decideNPolPlanes(casacore::Bool checkwithMS);
988 :
989 : //returns if mLocation_p is set (= non-default casacore::MPosition)
990 : casacore::Bool nonDefaultLocation();
991 :
992 : // implementation of mapExtent
993 : casacore::Bool getMapExtent(const casacore::MeasurementSet &ms,
994 : const casacore::String &referenceFrame, const casacore::String &movingSource,
995 : const casacore::String &pointingColumn, casacore::Vector<casacore::Double> ¢er, casacore::Vector<casacore::Double> &blc,
996 : casacore::Vector<casacore::Double> &trc, casacore::Vector<casacore::Double> &extent);
997 :
998 : ComponentList* componentList_p;
999 :
1000 : casacore::String scaleMethod_p; // "nscales" or "uservector"
1001 : casacore::Int nscales_p;
1002 : casacore::Int ntaylor_p;
1003 : casacore::Double reffreq_p;
1004 : casacore::Bool useNewMTFT_p;
1005 : casacore::Vector<casacore::Float> userScaleSizes_p;
1006 : casacore::Bool scaleInfoValid_p; // This means that we have set the information, not the scale beams
1007 : casacore::Float smallScaleBias_p; //ms-clean
1008 : casacore::Int nmodels_p;
1009 : // Everything here must be a real class since we make, handle and
1010 : // destroy these.
1011 : casacore::Block<casacore::CountedPtr<casacore::PagedImage<casacore::Float> > > images_p;
1012 : casacore::Block<casacore::CountedPtr<casacore::PagedImage<casacore::Float> > > masks_p;
1013 : casacore::Block<casacore::CountedPtr<casacore::PagedImage<casacore::Float> > > fluxMasks_p;
1014 : casacore::Block<casacore::CountedPtr<casacore::PagedImage<casacore::Float> > > residuals_p;
1015 :
1016 : // Freq frame is good and valid conversions can be done (or not)
1017 : casacore::Bool freqFrameValid_p;
1018 :
1019 : // Preferred complex polarization representation
1020 : StokesImageUtil::PolRep polRep_p;
1021 :
1022 : //Whether to use model column or use it in memory on the fly
1023 : casacore::Bool useModelCol_p;
1024 :
1025 : //Force single precision always
1026 : casacore::Bool singlePrec_p;
1027 : //sink used to store history mainly
1028 : casacore::LogSink logSink_p;
1029 :
1030 :
1031 : //
1032 : // Objects required for pointing correction (ftmachine=PBWProject)
1033 : //
1034 : EPJones *epJ;
1035 : casacore::String epJTableName_p, cfCacheDirName_p;
1036 : casacore::Bool doPointing, doPBCorr, psTermOn_p, aTermOn_p, mTermOn_p, wbAWP_p, conjBeams_p;
1037 : //SimplePlotterPtr plotter_p;
1038 : casacore::Record interactiveState_p;
1039 :
1040 : //Track moving source stuff
1041 : casacore::Bool doTrackSource_p;
1042 : casacore::MDirection trackDir_p;
1043 : casacore::String pointingDirCol_p;
1044 : VisImagingWeight imwgt_p;
1045 :
1046 : // viewer connection
1047 : int clean_panel_p;
1048 : int image_id_p;
1049 : int mask_id_p;
1050 : int prev_image_id_p;
1051 : int prev_mask_id_p;
1052 : //numthreads
1053 : casacore::Int numthreads_p;
1054 : casacore::Bool avoidTempLatt_p;
1055 : casacore::String projection_p;
1056 : };
1057 :
1058 :
1059 : } //# NAMESPACE CASA - END
1060 :
1061 : #endif
|