File: modelrenderer.cpp

package info (click to toggle)
wsclean 2.8-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye, sid
  • size: 2,196 kB
  • sloc: cpp: 34,504; ansic: 234; python: 174; makefile: 10
file content (312 lines) | stat: -rw-r--r-- 11,796 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312

#include <iostream>

#include "fitsreader.h"
#include "modelrenderer.h"
#include "model/model.h"
#include "units/imagecoordinates.h"
#include "uvector.h"
#include "fftconvolver.h"
#include "fftwmanager.h"

template<typename T>
T ModelRenderer::gaus(T x, T sigma)
{
	long double xi = x / sigma;
	return exp(T(-0.5) * xi * xi);// / (sigma * sqrt(T(2.0) * M_PIl));
}

/** Restore a circular beam*/
void ModelRenderer::Restore(double* imageData, size_t imageWidth, size_t imageHeight, const Model& model, long double beamSize, long double startFrequency, long double endFrequency, PolarizationEnum polarization)
{
	// Using the FWHM formula for a Gaussian:
	long double sigma = beamSize / (2.0L * sqrtl(2.0L * logl(2.0L)));
	
	int boundingBoxSize = ceil(sigma * 20.0 / std::min(_pixelScaleL, _pixelScaleM));
	for(Model::const_iterator src=model.begin(); src!=model.end(); ++src)
	{
		for(ModelSource::const_iterator comp=src->begin(); comp!=src->end(); ++comp)
		{
			long double
				posRA = comp->PosRA(),
				posDec = comp->PosDec(),
				sourceL, sourceM;
			ImageCoordinates::RaDecToLM(posRA, posDec, _phaseCentreRA, _phaseCentreDec, sourceL, sourceM);
			const SpectralEnergyDistribution &sed = comp->SED();
			const long double intFlux = sed.IntegratedFlux(startFrequency, endFrequency, polarization);
			
			//std::cout << "Source: " << comp->PosRA() << "," << comp->PosDec() << " Phase centre: " << _phaseCentreRA << "," << _phaseCentreDec << " beamsize: " << beamSize << "\n";
				
			int sourceX, sourceY;
			ImageCoordinates::LMToXY<long double>(sourceL-_phaseCentreDL, sourceM-_phaseCentreDM, _pixelScaleL, _pixelScaleM, imageWidth, imageHeight, sourceX, sourceY);
			//std::cout << "Adding source " << comp->Name() << " at " << sourceX << "," << sourceY << " of "
			//	<< intFlux << " Jy ("
			//	<< startFrequency/1000000.0 << "-" << endFrequency/1000000.0 << " MHz).\n";
			int
				xLeft = sourceX - boundingBoxSize,
				xRight = sourceX + boundingBoxSize,
				yTop = sourceY - boundingBoxSize,
				yBottom = sourceY + boundingBoxSize;
			if(xLeft < 0) xLeft = 0;
			if(xLeft > (int) imageWidth) xLeft = (int) imageWidth;
			if(xRight < xLeft) xRight = xLeft;
			if(xRight > (int) imageWidth) xRight = (int) imageWidth;
			if(yTop < 0) yTop = 0;
			if(yTop > (int) imageHeight) yTop = (int) imageHeight;
			if(yBottom < yTop) yBottom = yTop;
			if(yBottom > (int) imageHeight) yBottom = (int) imageHeight;
			
			for(int y=yTop; y!=yBottom; ++y)
			{
				double *imageDataPtr = imageData + y*imageWidth+xLeft;
				for(int x=xLeft; x!=xRight; ++x)
				{
					long double l, m;
					ImageCoordinates::XYToLM<long double>(x, y, _pixelScaleL, _pixelScaleM, imageWidth, imageHeight, l, m);
					l += _phaseCentreDL; m += _phaseCentreDM;
					long double dist = sqrt((l-sourceL)*(l-sourceL) + (m-sourceM)*(m-sourceM));
					long double g = gaus(dist, sigma);
					(*imageDataPtr) += double(g * intFlux);
					++imageDataPtr;
				}
			}
		}
	}
}

void ModelRenderer::renderGaussianComponent(double* imageData, size_t imageWidth, size_t imageHeight, long double posRA, long double posDec, long double gausMaj, long double gausMin, long double gausPA, long double flux)
{
	// Using the FWHM formula for a Gaussian:
	long double sigmaMaj = gausMaj / (2.0L * sqrtl(2.0L * logl(2.0L)));
	long double sigmaMin = gausMin / (2.0L * sqrtl(2.0L * logl(2.0L)));
	// TODO this won't work for non-equally spaced dimensions
	long double minPixelScale = std::min(_pixelScaleL, _pixelScaleM);
	// It is better to actually integrate the flux under the gaussian after
	// gridding, since a small gaussian might blow up with the following way.
	// double factor = 2.0L * M_PI * sigmaMaj * sigmaMin / (minPixelScale * minPixelScale);
	// flux /= factor;
	
	// Make rotation matrix
	long double transf[4];
	// Position angle is angle from North: 
	long double angle = gausPA+0.5*M_PI;
	transf[2] = sin(angle);
	transf[0] = cos(angle);
	transf[1] = -transf[2];
	transf[3] = transf[0];
	double sigmaMax = std::max(std::fabs(sigmaMaj * transf[0]), std::fabs(sigmaMaj * transf[1]));
	// Multiply with scaling matrix to make variance 1.
	transf[0] = transf[0] / sigmaMaj;
	transf[1] = transf[1] / sigmaMaj;
	transf[2] = transf[2] / sigmaMin;
	transf[3] = transf[3] / sigmaMin;
	int boundingBoxSize = ceil(sigmaMax * 20.0 / minPixelScale);
	long double sourceL, sourceM;
	ImageCoordinates::RaDecToLM(posRA, posDec, _phaseCentreRA, _phaseCentreDec, sourceL, sourceM);

	// Calculate the bounding box
	int sourceX, sourceY;
	ImageCoordinates::LMToXY<long double>(sourceL-_phaseCentreDL, sourceM-_phaseCentreDM, _pixelScaleL, _pixelScaleM, imageWidth, imageHeight, sourceX, sourceY);
	int
		xLeft = sourceX - boundingBoxSize,
		xRight = sourceX + boundingBoxSize,
		yTop = sourceY - boundingBoxSize,
		yBottom = sourceY + boundingBoxSize;
	if(xLeft < 0) xLeft = 0;
	if(xLeft > (int) imageWidth) xLeft = (int) imageWidth;
	if(xRight < xLeft) xRight = xLeft;
	if(xRight > (int) imageWidth) xRight = (int) imageWidth;
	if(yTop < 0) yTop = 0;
	if(yTop > (int) imageHeight) yTop = (int) imageHeight;
	if(yBottom < yTop) yBottom = yTop;
	if(yBottom > (int) imageHeight) yBottom = (int) imageHeight;
	
	ao::uvector<double> values;
	double fluxSum = 0.0;
	for(int y=yTop; y!=yBottom; ++y)
	{
		for(int x=xLeft; x!=xRight; ++x)
		{
			long double l, m;
			ImageCoordinates::XYToLM<long double>(x, y, _pixelScaleL, _pixelScaleM, imageWidth, imageHeight, l, m);
			l += _phaseCentreDL; m += _phaseCentreDM;
			long double
				lTransf = (l-sourceL)*transf[0] + (m-sourceM)*transf[1],
				mTransf = (l-sourceL)*transf[2] + (m-sourceM)*transf[3];
			long double dist = sqrt(lTransf*lTransf + mTransf*mTransf);
			long double v = gaus(dist, 1.0L) * flux;
			fluxSum += double(v);
			values.emplace_back(v);
		}
	}
	const double* iter = values.data();
	double factor = flux / fluxSum;
	for(int y=yTop; y!=yBottom; ++y)
	{
		double *imageDataPtr = imageData + y*imageWidth+xLeft;
		for(int x=xLeft; x!=xRight; ++x)
		{
			(*imageDataPtr) += *iter * factor;
			++imageDataPtr;
			++iter;
		}
	}
}

void ModelRenderer::renderPointComponent(double* imageData, size_t imageWidth, size_t imageHeight, long double posRA, long double posDec, long double flux)
{
	long double sourceL, sourceM;
	ImageCoordinates::RaDecToLM(posRA, posDec, _phaseCentreRA, _phaseCentreDec, sourceL, sourceM);
	sourceL -= _phaseCentreDL; sourceM -= _phaseCentreDM;
	
	int sourceX, sourceY;
	ImageCoordinates::LMToXY<long double>(sourceL, sourceM, _pixelScaleL, _pixelScaleM, imageWidth, imageHeight, sourceX, sourceY);
	
	if(sourceX >= 0 && sourceX < (int) imageWidth && sourceY >= 0 && sourceY < (int) imageHeight)
	{
		double *imageDataPtr = imageData + sourceY*imageWidth + sourceX;
		(*imageDataPtr) += double(flux);
	}
}

/** Restore a model with an elliptical beam */
void ModelRenderer::Restore(double* imageData, size_t imageWidth, size_t imageHeight, const Model& model, long double beamMaj, long double beamMin, long double beamPA, long double startFrequency, long double endFrequency, PolarizationEnum polarization)
{
	ao::uvector<double> renderedWithoutBeam(imageWidth * imageHeight, 0.0);
	renderModel(renderedWithoutBeam.data(), imageWidth, imageHeight, model, startFrequency, endFrequency, polarization);
	Restore(imageData, renderedWithoutBeam.data(), imageWidth, imageHeight, beamMaj, beamMin, beamPA, _pixelScaleL, _pixelScaleM);
}

/**
 * Restore a diffuse image (e.g. produced with multi-scale clean)
 */
void ModelRenderer::Restore(double* imageData, const double* modelData, size_t imageWidth, size_t imageHeight, long double beamMaj, long double beamMin, long double beamPA, long double pixelScaleL, long double pixelScaleM)
{
	if(beamMaj == 0.0 && beamMin == 0.0)
	{
		for(size_t j=0; j!=imageWidth*imageHeight; ++j)
			imageData[j] += modelData[j];
	}
	else {
		// Using the FWHM formula for a Gaussian:
		long double sigmaMaj = beamMaj / (2.0L * sqrtl(2.0L * logl(2.0L)));
		long double sigmaMin = beamMin / (2.0L * sqrtl(2.0L * logl(2.0L)));
		
		// Make rotation matrix
		long double transf[4];
		// Position angle is angle from North: 
		long double angle = beamPA+0.5*M_PI;
		transf[2] = sin(angle);
		transf[0] = cos(angle);
		transf[1] = -transf[2];
		transf[3] = transf[0];
		double sigmaMax = std::max(std::fabs(sigmaMaj * transf[0]), std::fabs(sigmaMaj * transf[1]));
		// Multiple with scaling matrix to make variance 1.
		transf[0] = transf[0] / sigmaMaj;
		transf[1] = transf[1] / sigmaMaj;
		transf[2] = transf[2] / sigmaMin;
		transf[3] = transf[3] / sigmaMin;
		
		size_t minDimension = std::min(imageWidth, imageHeight);
		size_t boundingBoxSize = std::min<size_t>(ceil(sigmaMax * 40.0 / std::min(pixelScaleL, pixelScaleM)), minDimension);
		if(boundingBoxSize%2!=0) ++boundingBoxSize;
		ao::uvector<double> kernel(boundingBoxSize*boundingBoxSize);
		typename ao::uvector<double>::iterator i=kernel.begin();
		for(size_t y=0; y!=boundingBoxSize; ++y)
		{
			for(size_t x=0; x!=boundingBoxSize; ++x)
			{
				long double l, m;
				ImageCoordinates::XYToLM<long double>(x, y, pixelScaleL, pixelScaleM, boundingBoxSize, boundingBoxSize, l, m);
				long double
					lTransf = l*transf[0] + m*transf[1],
					mTransf = l*transf[2] + m*transf[3];
				long double dist = sqrt(lTransf*lTransf + mTransf*mTransf);
				*i = gaus(dist, (long double) 1.0);
				++i;
			}
		}
		
		ao::uvector<double> convolvedModel(modelData, modelData+imageWidth*imageHeight);
		
		FFTWManager fftw;
		FFTConvolver::Convolve(fftw, convolvedModel.data(), imageWidth, imageHeight, kernel.data(), boundingBoxSize);
		for(size_t j=0; j!=imageWidth*imageHeight; ++j)
			imageData[j] += convolvedModel[j];
	}
}

/**
* Render without beam convolution, such that each point-source is one pixel.
*/
void ModelRenderer::renderModel(double* imageData, size_t imageWidth, size_t imageHeight, const Model& model, long double startFrequency, long double endFrequency, PolarizationEnum polarization)
{
	for(Model::const_iterator src=model.begin(); src!=model.end(); ++src)
	{
		for(ModelSource::const_iterator comp=src->begin(); comp!=src->end(); ++comp)
		{
			const long double
				posRA = comp->PosRA(),
				posDec = comp->PosDec(),
				intFlux = comp->SED().IntegratedFlux(startFrequency, endFrequency, polarization);
			
			if(comp->Type() == ModelComponent::GaussianSource)
			{
				long double
					gausMaj = comp->MajorAxis(),
					gausMin = comp->MinorAxis(),
					gausPA = comp->PositionAngle();
				renderGaussianComponent(imageData, imageWidth, imageHeight, posRA, posDec, gausMaj, gausMin, gausPA, intFlux);
			}
			else
				renderPointComponent(imageData, imageWidth, imageHeight, posRA, posDec, intFlux);
		}
	}
}

void ModelRenderer::RenderInterpolatedSource(double* image, size_t width, size_t height, double flux, double x, double y)
{
	ao::uvector<double>
		hSinc(std::min<size_t>(width,128)+1),
		vSinc(std::min<size_t>(height,128)+1);

	int midH = hSinc.size()/2;
	int midV = vSinc.size()/2;
	
	double xr = x - floor(x);
	double yr = y - floor(y);
	
	for(size_t i=0; i!=hSinc.size(); ++i)
	{
		double xi = (int(i) - midH - xr) * M_PI;
		if(xi == 0.0)
			hSinc[i] = 1.0;
		else
			hSinc[i] = sin(xi) / xi;
	}
	for(size_t i=0; i!=vSinc.size(); ++i)
	{
		double yi = (int(i) - midV - yr) * M_PI;
		if(yi == 0.0)
			vSinc[i] = 1.0;
		else
			vSinc[i] = sin(yi) / yi;
	}
	
	size_t
		xOffset = floor(x) - midH,
		yOffset = floor(y) - midV,
		startX = std::max<int>(xOffset, 0),
		startY = std::max<int>(yOffset, 0),
		endX = std::min<size_t>(xOffset + hSinc.size(), width),
		endY = std::min<size_t>(yOffset + vSinc.size(), height);
	for(size_t yi=startY; yi!=endY; ++yi)
	{
		double* ptr = &image[yi * width];
		double vFlux = flux * vSinc[yi - yOffset];
		for(size_t xi=startX; xi!=endX; ++xi)
		{
			ptr[xi] += vFlux * hSinc[xi - xOffset];
		}
	}
}