A
lgorithmique
N
umérique
D
istribuée
Public GIT Repository
projects
/
simgrid.git
/ blobdiff
commit
grep
author
committer
pickaxe
?
search:
re
summary
|
shortlog
|
log
|
commit
|
commitdiff
|
tree
raw
|
inline
| side by side
Update copyright lines with new year.
[simgrid.git]
/
src
/
smpi
/
colls
/
bcast
/
bcast-scatter-LR-allgather.cpp
diff --git
a/src/smpi/colls/bcast/bcast-scatter-LR-allgather.cpp
b/src/smpi/colls/bcast/bcast-scatter-LR-allgather.cpp
index
66ac6bb
..
7148af8
100644
(file)
--- a/
src/smpi/colls/bcast/bcast-scatter-LR-allgather.cpp
+++ b/
src/smpi/colls/bcast/bcast-scatter-LR-allgather.cpp
@@
-1,10
+1,10
@@
-/* Copyright (c) 2013-20
17
. The SimGrid Team. All rights reserved. */
+/* Copyright (c) 2013-20
20
. The SimGrid Team. All rights reserved. */
/* This program is free software; you can redistribute it and/or modify it
* under the terms of the license (GNU LGPL) which comes with this package. */
/* This program is free software; you can redistribute it and/or modify it
* under the terms of the license (GNU LGPL) which comes with this package. */
-#include "../colls_private.h"
-#include "s
rc/smpi/s
mpi_status.hpp"
+#include "../colls_private.h
pp
"
+#include "smpi_status.hpp"
/*****************************************************************************
/*****************************************************************************
@@
-64,21
+64,20
@@
SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
* Descrp: broadcasts using a scatter followed by LR allgather.
* Descrp: broadcasts using a scatter followed by LR allgather.
- * Auth
e
r: MPIH / modified by Ahmad Faraj
+ * Auth
o
r: MPIH / modified by Ahmad Faraj
****************************************************************************/
****************************************************************************/
-namespace simgrid{
-namespace smpi{
-int
-Coll_bcast_scatter_LR_allgather::bcast(void *buff, int count,
- MPI_Datatype data_type, int root,
- MPI_Comm comm)
+namespace simgrid {
+namespace smpi {
+int bcast__scatter_LR_allgather(void *buff, int count,
+ MPI_Datatype data_type, int root,
+ MPI_Comm comm)
{
MPI_Aint extent;
MPI_Status status;
int i, src, dst, rank, num_procs;
int mask, relative_rank, curr_size, recv_size, send_size, nbytes;
{
MPI_Aint extent;
MPI_Status status;
int i, src, dst, rank, num_procs;
int mask, relative_rank, curr_size, recv_size, send_size, nbytes;
- int scatter_size, left, right, next_src
, *recv_counts, *disps
;
+ int scatter_size, left, right, next_src;
int tag = COLL_TAG_BCAST;
rank = comm->rank();
int tag = COLL_TAG_BCAST;
rank = comm->rank();
@@
-87,7
+86,7
@@
Coll_bcast_scatter_LR_allgather::bcast(void *buff, int count,
nbytes = extent * count;
nbytes = extent * count;
- scatter_size = (nbytes + num_procs - 1) / num_procs; // ceiling division
+ scatter_size = (nbytes + num_procs - 1) / num_procs; // ceiling division
curr_size = (rank == root) ? nbytes : 0; // root starts with all the data
relative_rank = (rank >= root) ? rank - root : rank - root + num_procs;
curr_size = (rank == root) ? nbytes : 0; // root starts with all the data
relative_rank = (rank >= root) ? rank - root : rank - root + num_procs;
@@
-103,7
+102,7
@@
Coll_bcast_scatter_LR_allgather::bcast(void *buff, int count,
// allows you to post a larger recv.
if (recv_size <= 0)
curr_size = 0; // this process doesn't receive any data
// allows you to post a larger recv.
if (recv_size <= 0)
curr_size = 0; // this process doesn't receive any data
- // because of uneven division
+ // because of uneven division
else {
Request::recv((char *) buff + relative_rank * scatter_size, recv_size,
MPI_BYTE, src, tag, comm, &status);
else {
Request::recv((char *) buff + relative_rank * scatter_size, recv_size,
MPI_BYTE, src, tag, comm, &status);
@@
-123,7
+122,7
@@
Coll_bcast_scatter_LR_allgather::bcast(void *buff, int count,
while (mask > 0) {
if (relative_rank + mask < num_procs) {
send_size = curr_size - scatter_size * mask;
while (mask > 0) {
if (relative_rank + mask < num_procs) {
send_size = curr_size - scatter_size * mask;
- // mask is also the size of this process's subtree
+ // mask is also the size of this process's subtree
if (send_size > 0) {
dst = rank + mask;
if (send_size > 0) {
dst = rank + mask;
@@
-139,8
+138,8
@@
Coll_bcast_scatter_LR_allgather::bcast(void *buff, int count,
}
// done scatter now do allgather
}
// done scatter now do allgather
-
recv_counts = (int *) xbt_malloc(sizeof(int) * num_procs)
;
-
disps = (int *) xbt_malloc(sizeof(int) * num_procs)
;
+
int* recv_counts = new int[num_procs]
;
+
int* disps = new int[num_procs]
;
for (i = 0; i < num_procs; i++) {
recv_counts[i] = nbytes - i * scatter_size;
for (i = 0; i < num_procs; i++) {
recv_counts[i] = nbytes - i * scatter_size;
@@
-172,9
+171,8
@@
Coll_bcast_scatter_LR_allgather::bcast(void *buff, int count,
next_src = (num_procs + next_src - 1) % num_procs;
}
next_src = (num_procs + next_src - 1) % num_procs;
}
-
- free(recv_counts);
- free(disps);
+ delete[] recv_counts;
+ delete[] disps;
return MPI_SUCCESS;
}
return MPI_SUCCESS;
}