Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
H
HPC_Imaging
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Package registry
Container registry
Model registry
Operate
Environments
Terraform modules
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Claudio Gheller
HPC_Imaging
Commits
a22b9970
Commit
a22b9970
authored
3 years ago
by
Luca Tornatore
Browse files
Options
Downloads
Patches
Plain Diff
some sparse changes and compilation fixes
parent
d05e7f1e
No related branches found
No related tags found
No related merge requests found
Changes
3
Show whitespace changes
Inline
Side-by-side
Showing
3 changed files
gridding.c
+83
-73
83 additions, 73 deletions
gridding.c
init.c
+156
-142
156 additions, 142 deletions
init.c
w-stacking.cu
+2
-2
2 additions, 2 deletions
w-stacking.cu
with
241 additions
and
217 deletions
gridding.c
+
83
−
73
View file @
a22b9970
...
@@ -36,7 +36,8 @@ void initialize_array(){
...
@@ -36,7 +36,8 @@ void initialize_array(){
histo_send
=
(
long
*
)
calloc
(
nsectors
+
1
,
sizeof
(
long
));
histo_send
=
(
long
*
)
calloc
(
nsectors
+
1
,
sizeof
(
long
));
int
*
boundary
=
(
int
*
)
calloc
(
metaData
.
Nmeasures
,
sizeof
(
int
));
int
*
boundary
=
(
int
*
)
calloc
(
metaData
.
Nmeasures
,
sizeof
(
int
));
double
uuh
,
vvh
;
double
vvh
;
for
(
long
iphi
=
0
;
iphi
<
metaData
.
Nmeasures
;
iphi
++
)
for
(
long
iphi
=
0
;
iphi
<
metaData
.
Nmeasures
;
iphi
++
)
{
{
boundary
[
iphi
]
=
-
1
;
boundary
[
iphi
]
=
-
1
;
...
@@ -47,8 +48,10 @@ void initialize_array(){
...
@@ -47,8 +48,10 @@ void initialize_array(){
double
downdist
=
vvh
-
(
double
)(
binphi
*
yaxis
)
*
dx
;
double
downdist
=
vvh
-
(
double
)(
binphi
*
yaxis
)
*
dx
;
//
//
histo_send
[
binphi
]
++
;
histo_send
[
binphi
]
++
;
if
(
updist
<
w_supporth
&&
updist
>=
0
.
0
)
{
histo_send
[
binphi
+
1
]
++
;
boundary
[
iphi
]
=
binphi
+
1
;};
if
(
updist
<
w_supporth
&&
updist
>=
0
.
0
)
if
(
downdist
<
w_supporth
&&
binphi
>
0
&&
downdist
>=
0
.
0
)
{
histo_send
[
binphi
-
1
]
++
;
boundary
[
iphi
]
=
binphi
-
1
;};
{
histo_send
[
binphi
+
1
]
++
;
boundary
[
iphi
]
=
binphi
+
1
;};
if
(
downdist
<
w_supporth
&&
binphi
>
0
&&
downdist
>=
0
.
0
)
{
histo_send
[
binphi
-
1
]
++
;
boundary
[
iphi
]
=
binphi
-
1
;};
}
}
sectorarray
=
(
long
**
)
malloc
((
nsectors
+
1
)
*
sizeof
(
long
*
));
sectorarray
=
(
long
**
)
malloc
((
nsectors
+
1
)
*
sizeof
(
long
*
));
...
@@ -66,8 +69,10 @@ void initialize_array(){
...
@@ -66,8 +69,10 @@ void initialize_array(){
double
downdist
=
vvh
-
(
double
)(
binphi
*
yaxis
)
*
dx
;
double
downdist
=
vvh
-
(
double
)(
binphi
*
yaxis
)
*
dx
;
sectorarray
[
binphi
][
counter
[
binphi
]]
=
iphi
;
sectorarray
[
binphi
][
counter
[
binphi
]]
=
iphi
;
counter
[
binphi
]
++
;
counter
[
binphi
]
++
;
if
(
updist
<
w_supporth
&&
updist
>=
0
.
0
)
{
sectorarray
[
binphi
+
1
][
counter
[
binphi
+
1
]]
=
iphi
;
counter
[
binphi
+
1
]
++
;};
if
(
updist
<
w_supporth
&&
updist
>=
0
.
0
)
if
(
downdist
<
w_supporth
&&
binphi
>
0
&&
downdist
>=
0
.
0
)
{
sectorarray
[
binphi
-
1
][
counter
[
binphi
-
1
]]
=
iphi
;
counter
[
binphi
-
1
]
++
;};
{
sectorarray
[
binphi
+
1
][
counter
[
binphi
+
1
]]
=
iphi
;
counter
[
binphi
+
1
]
++
;};
if
(
downdist
<
w_supporth
&&
binphi
>
0
&&
downdist
>=
0
.
0
)
{
sectorarray
[
binphi
-
1
][
counter
[
binphi
-
1
]]
=
iphi
;
counter
[
binphi
-
1
]
++
;};
}
}
...
@@ -235,7 +240,7 @@ void gridding_data(){
...
@@ -235,7 +240,7 @@ void gridding_data(){
int
target_rank
=
(
int
)(
isector
%
size
);
int
target_rank
=
(
int
)(
isector
%
size
);
#ifdef ONE_SIDE
#ifdef ONE_SIDE
printf
(
"One Side communication active
\n
"
);
// MPI_Win_lock(MPI_LOCK_SHARED,target_rank,0,slabwin);
// MPI_Win_lock(MPI_LOCK_SHARED,target_rank,0,slabwin);
// MPI_Accumulate(gridss,size_of_grid,MPI_DOUBLE,target_rank,0,size_of_grid,MPI_DOUBLE,MPI_SUM,slabwin);
// MPI_Accumulate(gridss,size_of_grid,MPI_DOUBLE,target_rank,0,size_of_grid,MPI_DOUBLE,MPI_SUM,slabwin);
// MPI_Win_unlock(target_rank,slabwin);
// MPI_Win_unlock(target_rank,slabwin);
...
@@ -247,8 +252,9 @@ void gridding_data(){
...
@@ -247,8 +252,9 @@ void gridding_data(){
reduce
(
isector
,
target_rank
);
// here the reduce is performed within every host
reduce
(
isector
,
target_rank
);
// here the reduce is performed within every host
if
(
Me
.
Nhosts
>
1
)
// here thre reduce is performed among hosts
{
// here the reduce is performed among hosts
MPI_Barrier
(
MPI_COMM_WORLD
);
MPI_Barrier
(
MPI_COMM_WORLD
);
int
Im_in_the_new_communicator
=
MPI_UNDEFINED
;
int
Im_in_the_new_communicator
=
MPI_UNDEFINED
;
...
@@ -290,15 +296,15 @@ void gridding_data(){
...
@@ -290,15 +296,15 @@ void gridding_data(){
MPI_Comm_free
(
&
Sector_Comm
);
MPI_Comm_free
(
&
Sector_Comm
);
}
}
}
#else // relates to #ifdef ONE_SIDE
//MPI_Put(gridss,size_of_grid,MPI_DOUBLE,target_rank,0,size_of_grid,MPI_DOUBLE,slabwin);
#else
MPI_Reduce
(
gridss
,
grid
,
size_of_grid
,
MPI_DOUBLE
,
MPI_SUM
,
target_rank
,
MPI_COMM_WORLD
);
MPI_Reduce
(
gridss
,
grid
,
size_of_grid
,
MPI_DOUBLE
,
MPI_SUM
,
target_rank
,
MPI_COMM_WORLD
);
#endif //ONE_SIDE
#endif
//
closes #ifdef
ONE_SIDE
#endif //USE_MPI
#endif
//
closes
USE_MPI
clock_gettime
(
CLOCK_MONOTONIC
,
&
finishk
);
clock_gettime
(
CLOCK_MONOTONIC
,
&
finishk
);
endk
=
clock
();
endk
=
clock
();
...
@@ -430,6 +436,10 @@ void reduce( int sector, int target_rank )
...
@@ -430,6 +436,10 @@ void reduce( int sector, int target_rank )
int
local_rank
=
Me
.
Rank
[
myHOST
];
int
local_rank
=
Me
.
Rank
[
myHOST
];
if
(
Me
.
Ranks_to_host
[
target_rank
]
==
Me
.
myhost
)
if
(
Me
.
Ranks_to_host
[
target_rank
]
==
Me
.
myhost
)
// exchange rank 0 with target rank
// in this way the following log2 alogorithm,
// which reduces to rank 0, will work for
// every target rank
{
{
int
r
=
0
;
int
r
=
0
;
...
...
This diff is collapsed.
Click to expand it.
init.c
+
156
−
142
View file @
a22b9970
...
@@ -278,18 +278,20 @@ void readMetaData(char fileLocal[1000]) {
...
@@ -278,18 +278,20 @@ void readMetaData(char fileLocal[1000]) {
{
{
if
(
(
file
.
pFile
=
fopen
(
fileLocal
,
"r"
))
!=
NULL
)
if
(
(
file
.
pFile
=
fopen
(
fileLocal
,
"r"
))
!=
NULL
)
{
{
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
Nmeasures
);
int
items
=
0
;
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
Nvis
);
items
+=
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
Nmeasures
);
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
freq_per_chan
);
items
+=
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
Nvis
);
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
polarisations
);
items
+=
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
freq_per_chan
);
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
Ntimes
);
items
+=
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
polarisations
);
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
dt
);
items
+=
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
Ntimes
);
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
thours
);
items
+=
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
dt
);
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
baselines
);
items
+=
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
thours
);
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
uvmin
);
items
+=
fscanf
(
file
.
pFile
,
"%ld"
,
&
metaData
.
baselines
);
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
uvmax
);
items
+=
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
uvmin
);
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
wmin
);
items
+=
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
uvmax
);
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
wmax
);
items
+=
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
wmin
);
items
+=
fscanf
(
file
.
pFile
,
"%lf"
,
&
metaData
.
wmax
);
fclose
(
file
.
pFile
);
fclose
(
file
.
pFile
);
}
}
else
else
...
@@ -315,6 +317,7 @@ void metaData_calculation() {
...
@@ -315,6 +317,7 @@ void metaData_calculation() {
metaData
.
Nvis
=
metaData
.
Nmeasures
*
metaData
.
freq_per_chan
*
metaData
.
polarisations
;
metaData
.
Nvis
=
metaData
.
Nmeasures
*
metaData
.
freq_per_chan
*
metaData
.
polarisations
;
// calculate the coordinates of the center
// calculate the coordinates of the center
#warning why it it not used?
double
uvshift
=
metaData
.
uvmin
/
(
metaData
.
uvmax
-
metaData
.
uvmin
);
double
uvshift
=
metaData
.
uvmin
/
(
metaData
.
uvmax
-
metaData
.
uvmin
);
if
(
global_rank
==
0
)
if
(
global_rank
==
0
)
...
@@ -330,9 +333,7 @@ void metaData_calculation() {
...
@@ -330,9 +333,7 @@ void metaData_calculation() {
startrow
=
global_rank
*
nm_pe
;
startrow
=
global_rank
*
nm_pe
;
if
(
global_rank
==
size
-
1
)
nm_pe
=
nm_pe
+
remaining
;
if
(
global_rank
==
size
-
1
)
nm_pe
=
nm_pe
+
remaining
;
long
Nmeasures_tot
=
metaData
.
Nmeasures
;
metaData
.
Nmeasures
=
nm_pe
;
metaData
.
Nmeasures
=
nm_pe
;
long
Nvis_tot
=
metaData
.
Nvis
;
metaData
.
Nvis
=
metaData
.
Nmeasures
*
metaData
.
freq_per_chan
*
metaData
.
polarisations
;
metaData
.
Nvis
=
metaData
.
Nmeasures
*
metaData
.
freq_per_chan
*
metaData
.
polarisations
;
metaData
.
Nweights
=
metaData
.
Nmeasures
*
metaData
.
polarisations
;
metaData
.
Nweights
=
metaData
.
Nmeasures
*
metaData
.
polarisations
;
...
@@ -391,7 +392,9 @@ void readData() {
...
@@ -391,7 +392,9 @@ void readData() {
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
{
{
fseek
(
file
.
pFile
,
startrow
*
sizeof
(
double
),
SEEK_SET
);
fseek
(
file
.
pFile
,
startrow
*
sizeof
(
double
),
SEEK_SET
);
fread
(
data
.
uu
,
metaData
.
Nmeasures
*
sizeof
(
double
),
1
,
file
.
pFile
);
size_t
res
=
fread
(
data
.
uu
,
sizeof
(
double
),
metaData
.
Nmeasures
,
file
.
pFile
);
if
(
res
!=
metaData
.
Nmeasures
)
printf
(
"an error occurred while reading file %s
\n
"
,
filename
);
fclose
(
file
.
pFile
);
fclose
(
file
.
pFile
);
}
}
else
else
...
@@ -405,7 +408,10 @@ void readData() {
...
@@ -405,7 +408,10 @@ void readData() {
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
{
{
fseek
(
file
.
pFile
,
startrow
*
sizeof
(
double
),
SEEK_SET
);
fseek
(
file
.
pFile
,
startrow
*
sizeof
(
double
),
SEEK_SET
);
fread
(
data
.
vv
,
metaData
.
Nmeasures
*
sizeof
(
double
),
1
,
file
.
pFile
);
size_t
res
=
fread
(
data
.
vv
,
sizeof
(
double
),
metaData
.
Nmeasures
,
file
.
pFile
);
if
(
res
!=
metaData
.
Nmeasures
)
printf
(
"an error occurred while reading file %s
\n
"
,
filename
);
fclose
(
file
.
pFile
);
fclose
(
file
.
pFile
);
}
}
else
else
...
@@ -419,7 +425,9 @@ void readData() {
...
@@ -419,7 +425,9 @@ void readData() {
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
{
{
fseek
(
file
.
pFile
,
startrow
*
sizeof
(
double
),
SEEK_SET
);
fseek
(
file
.
pFile
,
startrow
*
sizeof
(
double
),
SEEK_SET
);
fread
(
data
.
ww
,
metaData
.
Nmeasures
*
sizeof
(
double
),
1
,
file
.
pFile
);
size_t
res
=
fread
(
data
.
ww
,
sizeof
(
double
),
metaData
.
Nmeasures
,
file
.
pFile
);
if
(
res
!=
metaData
.
Nmeasures
)
printf
(
"an error occurred while reading file %s
\n
"
,
filename
);
fclose
(
file
.
pFile
);
fclose
(
file
.
pFile
);
}
}
else
else
...
@@ -433,7 +441,9 @@ void readData() {
...
@@ -433,7 +441,9 @@ void readData() {
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
{
{
fseek
(
file
.
pFile
,
startrow
*
metaData
.
polarisations
*
sizeof
(
float
),
SEEK_SET
);
fseek
(
file
.
pFile
,
startrow
*
metaData
.
polarisations
*
sizeof
(
float
),
SEEK_SET
);
fread
(
data
.
weights
,(
metaData
.
Nweights
)
*
sizeof
(
float
),
1
,
file
.
pFile
);
size_t
res
=
fread
(
data
.
weights
,
sizeof
(
float
),
metaData
.
Nweights
,
file
.
pFile
);
if
(
res
!=
metaData
.
Nweights
)
printf
(
"an error occurred while reading file %s
\n
"
,
filename
);
fclose
(
file
.
pFile
);
fclose
(
file
.
pFile
);
}
}
else
else
...
@@ -447,7 +457,9 @@ void readData() {
...
@@ -447,7 +457,9 @@ void readData() {
if
((
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
if
((
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
{
{
fseek
(
file
.
pFile
,
startrow
*
metaData
.
freq_per_chan
*
metaData
.
polarisations
*
sizeof
(
float
),
SEEK_SET
);
fseek
(
file
.
pFile
,
startrow
*
metaData
.
freq_per_chan
*
metaData
.
polarisations
*
sizeof
(
float
),
SEEK_SET
);
fread
(
data
.
visreal
,
metaData
.
Nvis
*
sizeof
(
float
),
1
,
file
.
pFile
);
size_t
res
=
fread
(
data
.
visreal
,
sizeof
(
float
),
metaData
.
Nvis
,
file
.
pFile
);
if
(
res
!=
metaData
.
Nvis
)
printf
(
"an error occurred while reading file %s
\n
"
,
filename
);
fclose
(
file
.
pFile
);
fclose
(
file
.
pFile
);
}
}
else
else
...
@@ -461,7 +473,9 @@ void readData() {
...
@@ -461,7 +473,9 @@ void readData() {
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
if
(
(
file
.
pFile
=
fopen
(
filename
,
"rb"
))
!=
NULL
)
{
{
fseek
(
file
.
pFile
,
startrow
*
metaData
.
freq_per_chan
*
metaData
.
polarisations
*
sizeof
(
float
),
SEEK_SET
);
fseek
(
file
.
pFile
,
startrow
*
metaData
.
freq_per_chan
*
metaData
.
polarisations
*
sizeof
(
float
),
SEEK_SET
);
fread
(
data
.
visimg
,
metaData
.
Nvis
*
sizeof
(
float
),
1
,
file
.
pFile
);
size_t
res
=
fread
(
data
.
visimg
,
sizeof
(
float
),
metaData
.
Nvis
,
file
.
pFile
);
if
(
res
!=
metaData
.
Nvis
)
printf
(
"an error occurred while reading file %s
\n
"
,
filename
);
fclose
(
file
.
pFile
);
fclose
(
file
.
pFile
);
}
}
else
else
...
...
This diff is collapsed.
Click to expand it.
w-stacking.cu
+
2
−
2
View file @
a22b9970
...
@@ -140,7 +140,7 @@ void wstack(
...
@@ -140,7 +140,7 @@ void wstack(
int
num_threads
)
int
num_threads
)
{
{
long
i
;
long
i
;
long
index
;
//
long index;
long
visindex
;
long
visindex
;
// initialize the convolution kernel
// initialize the convolution kernel
...
@@ -242,7 +242,7 @@ void wstack(
...
@@ -242,7 +242,7 @@ void wstack(
visindex
=
i
*
freq_per_chan
*
polarizations
;
visindex
=
i
*
freq_per_chan
*
polarizations
;
double
sum
=
0.0
;
//
double sum = 0.0;
int
j
,
k
;
int
j
,
k
;
//if (i%1000 == 0)printf("%ld\n",i);
//if (i%1000 == 0)printf("%ld\n",i);
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment