-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathdata_processing-mpi.c
225 lines (200 loc) · 9.03 KB
/
data_processing-mpi.c
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
/*
* Copyright (c) 2024 Marzia Rivi
*
* This file is part of RadioLensfit.
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 2 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, see <http://www.gnu.org/licenses/>.
*/
// data_processing_mpi.c
#ifdef USE_MPI
#include <mpi.h>
#endif
#ifdef _OPENMP
#include <omp.h>
#endif
#include <stdio.h>
#include <math.h>
#include <string.h>
#include "datatype.h"
#include "default_params.h"
#include "utils.h"
#include "galaxy_fitting.h"
#include "galaxy_visibilities.h"
#include "distributions.h"
#include "evaluate_uv_grid.h"
#include "source_extraction.h"
#include "data_processing.h"
#ifdef __cplusplus
extern "C" {
#endif
#ifdef USE_MPI
void data_processing_mpi(bool re_fitting, unsigned int *bad_list, int nprocs, int rank, int nsources, double len, unsigned int num_coords,
FILE *pFile, likelihood_params *par, double *l, double *m, double *gflux, double *gscale, double *ge1, double *ge2, double *SNR_vis,
double *sum_w, complexd *visGal, complexd *visSkyMod, complexd *visData,
float *sigma2_vis, bool *flag, double *uu_metres, double *vv_metres, double *ww_metres, complexd *temp_facet_visData,
double *temp_facet_sigma2, double *temp_sum, double *com_time, double *fitting_time, int *bad)
{
MPI_Status stat;
unsigned int gal;
double l0,m0;
int k, nbad;
long int my_g, ind;
double flux, mu, R_mu[nprocs];
nbad = 0;
unsigned int g = 0; // source global index
while (g < nsources)
{
k = 0; // source local index
my_g = -1; // if sources are finished, current task has my_g = -1 and has no source to fit
// source extraction of nproc consecutive sources that will be fitted each one by a different task -----------------------
for(int src = 0; src < nprocs && g < nsources; src++)
{
if (re_fitting) gal = bad_list[g];
else gal = g;
l0 = l[gal]; m0 = m[gal];
flux = gflux[gal];
mu = scale_mean(flux);
#ifndef SCALELENGTH_ON
R_mu[src] = exp(mu);
#else
R_mu[src] = gscale[gal];
#endif
unsigned int facet = facet_size(R_mu[src],len);
unsigned long int size = (unsigned long int) facet*facet;
if (rank == src) // proc src will fit the current source:
{
my_g = gal;
for (int nRo=1; nRo<par->numr; nRo++)
(par->rprior)[nRo] = rfunc(mu,R_STD,par->ro[nRo]); // set log(prior) for scalelength of source gal
// extract visibilities, sigma2 and weights from my MS (already summed in the facet) for source gal
par->facet = facet;
source_extraction(rank,facet,par,par->data,par->sigma2,sum_w,l0, m0, flux, R_mu[src], 0., 0., visSkyMod, visData, visGal, sigma2_vis, flag, num_coords, uu_metres, vv_metres, ww_metres, len);
int n = 1;
while (n<nprocs)
{
// collect and reduce facet vis, sigma2 and sum_w of the current source from the other procs (for their MS contribution)
*com_time -= MPI_Wtime();
MPI_Recv(temp_sum,size,MPI_DOUBLE,MPI_ANY_SOURCE,src,MPI_COMM_WORLD,&stat);
MPI_Recv(temp_facet_sigma2,size,MPI_DOUBLE,MPI_ANY_SOURCE,nprocs+src,MPI_COMM_WORLD,&stat);
MPI_Recv(temp_facet_visData,2*size,MPI_DOUBLE,MPI_ANY_SOURCE,2*nprocs+src,MPI_COMM_WORLD,&stat);
*com_time += MPI_Wtime();
for (unsigned long int i = 0; i<size; i++) sum_w[i] += temp_sum[i];
for (unsigned long int i = 0; i<size; i++) (par->sigma2)[i] += temp_facet_sigma2[i];
for (unsigned long int i = 0; i<size; i++)
{
(par->data)[i].real += temp_facet_visData[i].real;
(par->data)[i].imag += temp_facet_visData[i].imag;
}
n++;
}
}
else
{
// extract visibilities, sigma2 and sum_w from my MS (already summed in the facet) for source gal and send them to proc src
source_extraction(rank,facet,par,temp_facet_visData, temp_facet_sigma2,temp_sum,l0, m0, flux, R_mu[src], 0., 0., visSkyMod, visData, visGal, sigma2_vis, flag, num_coords, uu_metres, vv_metres, ww_metres, len);
*com_time -= MPI_Wtime();
MPI_Send(temp_sum,size,MPI_DOUBLE,src,src,MPI_COMM_WORLD);
MPI_Send(temp_facet_sigma2,size,MPI_DOUBLE,src,nprocs+src,MPI_COMM_WORLD);
MPI_Send(temp_facet_visData,2*size,MPI_DOUBLE,src,2*nprocs+src,MPI_COMM_WORLD);
*com_time += MPI_Wtime();
}
g++;
k++;
}
if (my_g >= 0)
{
// average facet visibilities (already summed within the cell)
average_facets(par->facet*par->facet, par->data, par->sigma2, sum_w);
// compute facet coordinates their number
par->ncoords = evaluate_facet_coords(par->uu, par->vv, len, par->facet, sum_w);
}
// source fitting of this task --------------------------------------------------------------------------------------------------
double start_fitting = MPI_Wtime();
double mes_e1, mes_e2, maxL;
double var_e1, var_e2, oneDimvar;
if (my_g >= 0)
{
source_fitting(rank, par, &mes_e1, &mes_e2, &var_e1, &var_e2, &oneDimvar, &maxL);
printf("rank %d: n. %d flux = %f: measured e = %f , %f \n",rank,my_g,gflux[my_g],mes_e1,mes_e2);
}
*fitting_time += MPI_Wtime() - start_fitting;
// removal of the fitted sources visibilities from the MS data and sky model ----------------------------------------------------------
double res[6]; // shape fitting results to be sent to the other procs
ind = g - k; // source global index
k = 0; // source local index
for (int src = 0 ; src < nprocs && ind < nsources; src++)
{
if (re_fitting) gal = bad_list[ind]; //source global index
else gal = ind;
l0 = l[gal]; m0 = m[gal];
flux = gflux[gal];
if (rank == k)
{
res[0] = mes_e1; res[1] = mes_e2;
res[2] = var_e1; res[3] = var_e2;
res[4] = oneDimvar; res[5] = maxL;
}
//MPI_Barrier(MPI_COMM_WORLD);
// Bcast shape results from proc k to the others
*com_time -= MPI_Wtime();
MPI_Bcast(res,6,MPI_DOUBLE,k,MPI_COMM_WORLD);
*com_time += MPI_Wtime();
if (res[5] <= -1e+10 || res[2] < VAR || res[3] < VAR || res[4] < VAR) // bad measurement
{
if (re_fitting)
{
if (rank == 0) fprintf(pFile, "%f | %f | %f | %e | %f | %f | %e | %e | %f | %f | %f \n",flux,ge1[gal],res[0],sqrt(res[2]),ge2[gal],res[1],sqrt(res[3]),res[4],SNR_vis[gal],l0/(ARCS2RAD),m0/(ARCS2RAD));
res[0] = 0.; res[1] = 0.; // remove round source model from original data
}
else bad_list[nbad] = ind; // each proc store the index of bad measurements to be fit again at the end
nbad++;
}
else
{
#ifdef _OPENMP
#pragma omp parallel for
#endif
for (unsigned int ch = 0; ch < par->nchannels; ch++)
{
unsigned long int ch_vis = (unsigned long int) ch*num_coords;
// remove current round source model from the sky model
data_galaxy_visibilities((par->spec)[ch], (par->wavenumbers)[ch], par->band_factor, par->acc_time, 0., 0., R_mu[k],
flux, l0, m0, num_coords, uu_metres, vv_metres, ww_metres, &(visGal[ch_vis]));
for (unsigned long int i = ch_vis; i<ch_vis+num_coords; i++)
{
visSkyMod[i].real -= visGal[i].real;
visSkyMod[i].imag -= visGal[i].imag;
}
// remove current source model fit from original data
data_galaxy_visibilities((par->spec)[ch], (par->wavenumbers)[ch], par->band_factor, par->acc_time, res[0], res[1], R_mu[k],
flux, l0, m0, num_coords, uu_metres, vv_metres, ww_metres, &(visGal[ch_vis]));
for (unsigned long int i = ch_vis; i<ch_vis+num_coords; i++)
{
visData[i].real -= visGal[i].real;
visData[i].imag -= visGal[i].imag;
}
}
if (rank == 0) fprintf(pFile, "%f | %f | %f | %e | %f | %f | %e | %e | %f | %f | %f \n",flux,ge1[gal],res[0],sqrt(res[2]),ge2[gal],res[1],sqrt(res[3]),res[4],SNR_vis[gal],l0/(ARCS2RAD),m0/(ARCS2RAD));
}
k++;
ind++;
}
}
*bad = nbad;
return;
}
#endif
#ifdef __cplusplus
}
#endif