Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F122579516
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
Sun, Jul 20, 17:32
Size
3 KB
Mime Type
text/x-c
Expires
Tue, Jul 22, 17:32 (1 d, 23 h)
Engine
blob
Format
Raw Data
Handle
27507213
Attached To
rSCINTROPARALLEL Poisson code for introduction to parallelism
grid.cc
View Options
/* -------------------------------------------------------------------------- */
#include "grid.hh"
/* -------------------------------------------------------------------------- */
#include <algorithm>
#if defined(_DEBUG_MPI)
#include <iomanip>
#endif
/* -------------------------------------------------------------------------- */
/* -------------------------------------------------------------------------- */
Grid
::
Grid
(
int
m
,
int
n
)
:
m_m
(
m
),
m_n
(
n
),
m_storage
(
m
*
n
)
{
clear
();
}
/* -------------------------------------------------------------------------- */
Grid
::~
Grid
()
{
#if defined(_MPI)
MPI_Type_free
(
&
column_t
);
MPI_Type_free
(
&
line_t
);
#endif
}
/* -------------------------------------------------------------------------- */
void
Grid
::
clear
()
{
std
::
fill
(
m_storage
.
begin
(),
m_storage
.
end
(),
0.
);
}
/* -------------------------------------------------------------------------- */
void
Grid
::
resize
(
int
m
,
int
n
)
{
m_m
=
m
;
m_n
=
n
;
m_storage
.
resize
(
m
*
n
);
}
/* -------------------------------------------------------------------------- */
int
Grid
::
m
()
const
{
return
m_m
;
}
int
Grid
::
n
()
const
{
return
m_n
;
}
/* -------------------------------------------------------------------------- */
#if defined(_MPI)
void
Grid
::
initializeHaloCommunications
(
MPI_Comm
&
m_communicator
,
const
std
::
array
<
int
,
2
>
&
ghosts
)
{
std
::
array
<
int
,
2
>
coords
;
int
prank
;
MPI_Type_vector
(
m_m
-
ghosts
[
0
],
1
,
m_n
,
MPI_FLOAT
,
&
column_t
);
MPI_Type_vector
(
1
,
m_n
-
ghosts
[
1
],
1
,
MPI_FLOAT
,
&
line_t
);
MPI_Type_commit
(
&
column_t
);
MPI_Type_commit
(
&
line_t
);
MPI_Comm_rank
(
m_communicator
,
&
prank
);
MPI_Cart_coords
(
m_communicator
,
prank
,
coords
.
size
(),
coords
.
data
());
// determining the rank of the neighbors
int
left
,
right
,
top
,
bottom
;
MPI_Cart_shift
(
m_communicator
,
0
,
1
,
&
left
,
&
right
);
MPI_Cart_shift
(
m_communicator
,
1
,
1
,
&
top
,
&
bottom
);
#if defined(_DEBUG_MPI)
static
bool
first
=
true
;
if
(
first
)
{
int
psize
;
MPI_Comm_size
(
m_communicator
,
&
psize
);
for
(
int
i
=
0
;
i
<
psize
;
++
i
)
{
if
(
prank
==
i
)
{
std
::
cerr
<<
prank
<<
std
::
setw
(
7
)
<<
top
<<
"
\n
"
<<
prank
<<
std
::
setw
(
8
)
<<
"^
\n
"
<<
prank
<<
std
::
setw
(
3
)
<<
left
<<
" < . > "
<<
right
<<
"
\n
"
<<
prank
<<
std
::
setw
(
8
)
<<
"v
\n
"
<<
prank
<<
std
::
setw
(
7
)
<<
bottom
<<
"
\n\n
"
;
}
MPI_Barrier
(
m_communicator
);
}
first
=
false
;
}
#endif
auto
addr
=
[
&
](
auto
&&
i
,
auto
&&
j
)
{
return
&
(
this
->
operator
()(
i
,
j
));
};
MPI_Request
*
rrequest
=
m_requests
.
data
();
MPI_Request
*
srequest
=
m_requests
.
data
()
+
4
;
static
int
tag
=
1
;
/// preparing receives
{
auto
send
=
addr
(
m_m
-
2
,
1
);
auto
recv
=
addr
(
0
,
1
);
MPI_Send_init
(
send
,
1
,
line_t
,
bottom
,
tag
,
m_communicator
,
srequest
);
MPI_Recv_init
(
recv
,
1
,
line_t
,
top
,
tag
,
m_communicator
,
rrequest
);
++
rrequest
;
++
srequest
;
++
tag
;
}
{
auto
send
=
addr
(
1
,
m_n
-
2
);
auto
recv
=
addr
(
1
,
0
);
MPI_Send_init
(
send
,
1
,
column_t
,
right
,
tag
,
m_communicator
,
srequest
);
MPI_Recv_init
(
recv
,
1
,
column_t
,
left
,
tag
,
m_communicator
,
rrequest
);
++
rrequest
;
++
srequest
;
++
tag
;
}
{
auto
send
=
addr
(
1
,
1
);
auto
recv
=
addr
(
m_m
-
1
,
1
);
MPI_Send_init
(
send
,
1
,
line_t
,
top
,
tag
,
m_communicator
,
srequest
);
MPI_Recv_init
(
recv
,
1
,
line_t
,
bottom
,
tag
,
m_communicator
,
rrequest
);
++
rrequest
;
++
srequest
;
++
tag
;
}
{
auto
send
=
addr
(
1
,
1
);
auto
recv
=
addr
(
1
,
m_n
-
1
);
MPI_Send_init
(
send
,
1
,
column_t
,
left
,
tag
,
m_communicator
,
srequest
);
MPI_Recv_init
(
recv
,
1
,
column_t
,
right
,
tag
,
m_communicator
,
rrequest
);
++
rrequest
;
++
srequest
;
++
tag
;
}
}
#endif
Event Timeline
Log In to Comment