Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F67070557
grid.cc
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Wed, Jun 19, 18:06
Size
4 KB
Mime Type
text/x-c
Expires
Fri, Jun 21, 18:06 (2 d)
Engine
blob
Format
Raw Data
Handle
18334949
Attached To
rSCMINICLUSTER SCITAS Mini Cluster
grid.cc
View Options
#include "grid.hh"
#include <QDebug>
#include <mpi.h>
typedef
enum
{
_north
,
_south
,
_east
,
_west
,
_north_east
,
_north_west
,
_south_east
,
_south_west
}
comm_tag_t
;
Grid
::~
Grid
()
{
delete
[]
_grid
;
MPI_Type_free
(
&
_column_t
);
MPI_Type_free
(
&
_line_t
);
MPI_Type_free
(
&
_corner_t
);
MPI_Request
*
req
=
_requests
;
for
(
int
i
=
0
;
i
<
8
*
2
;
++
i
,
++
req
)
{
MPI_Request_free
(
req
);
}
delete
[]
_requests
;
}
void
Grid
::
initalizeCommunicationRequests
(
MPI_Comm
&
comm
,
int
coords
[
2
])
{
_nb_requests
=
8
;
_requests
=
new
MPI_Request
[
_nb_requests
*
2
];
MPI_Type_vector
(
_dims
[
0
],
_ghosts
,
_dims
[
1
]
+
2
*
_ghosts
,
MPI_CHAR
,
&
_column_t
);
MPI_Type_vector
(
_ghosts
,
_dims
[
1
],
2
*
_ghosts
,
MPI_CHAR
,
&
_line_t
);
MPI_Type_vector
(
_ghosts
,
_ghosts
,
_dims
[
1
]
+
2
*
_ghosts
,
MPI_CHAR
,
&
_corner_t
);
MPI_Type_commit
(
&
_column_t
);
MPI_Type_commit
(
&
_line_t
);
MPI_Type_commit
(
&
_corner_t
);
int
prank
;
MPI_Comm_rank
(
comm
,
&
prank
);
auto
data_addr
=
[
this
](
int
i
,
int
j
)
->
char
*
{
if
(
i
<
0
)
{
i
=
_dims
[
0
]
+
2
*
_ghosts
+
i
;
}
if
(
j
<
0
)
{
j
=
_dims
[
1
]
+
2
*
_ghosts
+
j
;
}
return
_grid
+
i
*
(
_dims
[
1
]
+
2
*
_ghosts
)
+
j
;
};
MPI_Request
*
srequests
=
_requests
;
MPI_Request
*
rrequests
=
_requests
+
_nb_requests
;
{
/* south east -> north west */
int
nrank
;
int
rcoords
[
2
]
=
{
coords
[
0
]
-
1
,
coords
[
1
]
-
1
};
MPI_Cart_rank
(
comm
,
rcoords
,
&
nrank
);
char
*
send
=
data_addr
(
_ghosts
,
_ghosts
);
char
*
recv
=
data_addr
(
0
,
0
);
MPI_Send_init
(
send
,
1
,
_corner_t
,
nrank
,
_south_east
,
comm
,
srequests
);
MPI_Recv_init
(
recv
,
1
,
_corner_t
,
nrank
,
_north_west
,
comm
,
rrequests
);
++
srequests
;
++
rrequests
;
}
{
/* south west -> north east */
int
nrank
;
int
rcoords
[
2
]
=
{
coords
[
0
]
-
1
,
coords
[
1
]
+
1
};
MPI_Cart_rank
(
comm
,
rcoords
,
&
nrank
);
char
*
send
=
data_addr
(
_ghosts
,
-
2
*
_ghosts
);
char
*
recv
=
data_addr
(
0
,
-
_ghosts
);
MPI_Send_init
(
send
,
1
,
_corner_t
,
nrank
,
_south_west
,
comm
,
srequests
);
MPI_Recv_init
(
recv
,
1
,
_corner_t
,
nrank
,
_north_east
,
comm
,
rrequests
);
++
srequests
;
++
rrequests
;
}
{
/* north west -> south east */
int
nrank
;
int
rcoords
[
2
]
=
{
coords
[
0
]
+
1
,
coords
[
1
]
+
1
};
MPI_Cart_rank
(
comm
,
rcoords
,
&
nrank
);
char
*
send
=
data_addr
(
-
2
*
_ghosts
,
-
2
*
_ghosts
);
char
*
recv
=
data_addr
(
-
_ghosts
,
-
_ghosts
);
MPI_Send_init
(
send
,
1
,
_corner_t
,
nrank
,
_north_west
,
comm
,
srequests
);
MPI_Recv_init
(
recv
,
1
,
_corner_t
,
nrank
,
_south_east
,
comm
,
rrequests
);
++
srequests
;
++
rrequests
;
}
{
/* north_east -> south_west */
int
nrank
;
int
rcoords
[
2
]
=
{
coords
[
0
]
+
1
,
coords
[
1
]
-
1
};
MPI_Cart_rank
(
comm
,
rcoords
,
&
nrank
);
char
*
send
=
data_addr
(
-
2
*
_ghosts
,
_ghosts
);
char
*
recv
=
data_addr
(
-
_ghosts
,
0
);
MPI_Send_init
(
send
,
1
,
_corner_t
,
nrank
,
_north_east
,
comm
,
srequests
);
MPI_Recv_init
(
recv
,
1
,
_corner_t
,
nrank
,
_south_west
,
comm
,
rrequests
);
++
srequests
;
++
rrequests
;
}
// Sides
{
/* south -> north */
int
nrank
;
int
rcoords
[
2
]
=
{
coords
[
0
]
-
1
,
coords
[
1
]};
MPI_Cart_rank
(
comm
,
rcoords
,
&
nrank
);
char
*
send
=
data_addr
(
_ghosts
,
_ghosts
);
char
*
recv
=
data_addr
(
0
,
_ghosts
);
MPI_Send_init
(
send
,
1
,
_line_t
,
nrank
,
_north
,
comm
,
srequests
);
MPI_Recv_init
(
recv
,
1
,
_line_t
,
nrank
,
_south
,
comm
,
rrequests
);
++
srequests
;
++
rrequests
;
}
{
/* west -> east */
int
nrank
;
int
rcoords
[
2
]
=
{
coords
[
0
],
coords
[
1
]
+
1
};
MPI_Cart_rank
(
comm
,
rcoords
,
&
nrank
);
char
*
send
=
data_addr
(
_ghosts
,
-
2
*
_ghosts
);
char
*
recv
=
data_addr
(
_ghosts
,
-
_ghosts
);
MPI_Send_init
(
send
,
1
,
_column_t
,
nrank
,
_west
,
comm
,
srequests
);
MPI_Recv_init
(
recv
,
1
,
_column_t
,
nrank
,
_east
,
comm
,
rrequests
);
++
srequests
;
++
rrequests
;
}
{
/* north -> south */
int
nrank
;
int
rcoords
[
2
]
=
{
coords
[
0
]
+
1
,
coords
[
1
]};
MPI_Cart_rank
(
comm
,
rcoords
,
&
nrank
);
char
*
send
=
data_addr
(
-
2
*
_ghosts
,
_ghosts
);
char
*
recv
=
data_addr
(
-
_ghosts
,
_ghosts
);
MPI_Send_init
(
send
,
1
,
_line_t
,
nrank
,
_south
,
comm
,
srequests
);
MPI_Recv_init
(
recv
,
1
,
_line_t
,
nrank
,
_north
,
comm
,
rrequests
);
++
srequests
;
++
rrequests
;
}
{
/* east -> west */
int
nrank
;
int
rcoords
[
2
]
=
{
coords
[
0
],
coords
[
1
]
-
1
};
MPI_Cart_rank
(
comm
,
rcoords
,
&
nrank
);
char
*
send
=
data_addr
(
_ghosts
,
_ghosts
);
char
*
recv
=
data_addr
(
_ghosts
,
0
);
MPI_Send_init
(
send
,
1
,
_column_t
,
nrank
,
_east
,
comm
,
srequests
);
MPI_Recv_init
(
recv
,
1
,
_column_t
,
nrank
,
_west
,
comm
,
rrequests
);
++
srequests
;
++
rrequests
;
}
}
Event Timeline
Log In to Comment