Skip to content
GitLab
Menu
Projects
Groups
Snippets
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Sign in
Toggle navigation
Menu
Open sidebar
Podhorszki, Norbert
ADIOS2
Commits
0fb306c9
Commit
0fb306c9
authored
Mar 08, 2020
by
Ruonan Wang
Browse files
xgc-coupler-gene 3 way coupling test passed
parent
5c57e9f4
Changes
4
Hide whitespace changes
Inline
Side-by-side
source/adios2/engine/ssc/SscReader.cpp
View file @
0fb306c9
...
...
@@ -197,7 +197,6 @@ void SscReader::SyncMpiPattern()
m_MaxFilenameLength
,
m_RendezvousAppCount
,
CommAsMPI
(
m_Comm
));
m_MpiAllWritersGroup
=
m_MpiHandshake
.
GetAllWritersGroup
(
m_Name
);
m_StreamComm
=
m_MpiHandshake
.
GetStreamComm
(
m_Name
);
}
...
...
source/adios2/helper/adiosMpiHandshake.cpp
View file @
0fb306c9
...
...
@@ -278,7 +278,8 @@ MPI_Comm MpiHandshake::GetStreamComm(const std::string &filename)
MPI_Comm_group
(
MPI_COMM_WORLD
,
&
worldGroup
);
std
::
sort
(
allStreamRanks
.
begin
(),
allStreamRanks
.
end
());
MPI_Group
allWorkersGroup
;
MPI_Group_incl
(
worldGroup
,
allStreamRanks
.
size
(),
allStreamRanks
.
data
(),
&
allWorkersGroup
);
MPI_Group_incl
(
worldGroup
,
allStreamRanks
.
size
(),
allStreamRanks
.
data
(),
&
allWorkersGroup
);
MPI_Comm
streamComm
;
MPI_Comm_create_group
(
MPI_COMM_WORLD
,
allWorkersGroup
,
0
,
&
streamComm
);
return
streamComm
;
...
...
@@ -299,7 +300,8 @@ MPI_Group MpiHandshake::GetAllReadersGroup(const std::string &filename)
MPI_Group
worldGroup
;
MPI_Group
allReadersGroup
;
MPI_Comm_group
(
MPI_COMM_WORLD
,
&
worldGroup
);
MPI_Group_incl
(
worldGroup
,
allReaderRanks
.
size
(),
allReaderRanks
.
data
(),
&
allReadersGroup
);
MPI_Group_incl
(
worldGroup
,
allReaderRanks
.
size
(),
allReaderRanks
.
data
(),
&
allReadersGroup
);
return
allReadersGroup
;
}
...
...
@@ -318,7 +320,8 @@ MPI_Group MpiHandshake::GetAllWritersGroup(const std::string &filename)
MPI_Group
worldGroup
;
MPI_Group
allWritersGroup
;
MPI_Comm_group
(
MPI_COMM_WORLD
,
&
worldGroup
);
MPI_Group_incl
(
worldGroup
,
allWriterRanks
.
size
(),
allWriterRanks
.
data
(),
&
allWritersGroup
);
MPI_Group_incl
(
worldGroup
,
allWriterRanks
.
size
(),
allWriterRanks
.
data
(),
&
allWritersGroup
);
return
allWritersGroup
;
}
...
...
source/adios2/helper/adiosMpiHandshake.h
View file @
0fb306c9
...
...
@@ -96,7 +96,8 @@ public:
GetReaderMap
(
const
std
::
string
&
filename
);
/**
* Get the MPI communicator for the stream, containing all writer ranks and reader ranks that work on the stream filename
* Get the MPI communicator for the stream, containing all writer ranks and
* reader ranks that work on the stream filename
*
* @param filename: name of the staging stream
*
...
...
testing/adios2/engine/ssc/TestSscXgc3Way.cpp
View file @
0fb306c9
...
...
@@ -62,30 +62,15 @@ void coupler(const Dims &shape, const Dims &start, const Dims &count,
auto
c_to_g_var
=
c_to_g_io
.
DefineVariable
<
float
>
(
"c_to_g"
,
shape
,
start
,
count
);
std
::
cout
<<
" ======================== coupler before handshake "
<<
std
::
endl
;
adios2
::
Engine
x_to_c_engine
=
x_to_c_io
.
Open
(
"x_to_c"
,
adios2
::
Mode
::
Read
);
std
::
cout
<<
" ======================== "
<<
std
::
endl
;
std
::
cout
<<
" x_to_c handshake finished on coupler rank "
<<
worldRank
<<
std
::
endl
;
adios2
::
Engine
c_to_x_engine
=
c_to_x_io
.
Open
(
"c_to_x"
,
adios2
::
Mode
::
Write
);
std
::
cout
<<
" ======================== "
<<
std
::
endl
;
std
::
cout
<<
" c_to_x handshake finished on coupler rank "
<<
worldRank
<<
std
::
endl
;
adios2
::
Engine
c_to_g_engine
=
c_to_g_io
.
Open
(
"c_to_g"
,
adios2
::
Mode
::
Write
);
std
::
cout
<<
" ======================== "
<<
std
::
endl
;
std
::
cout
<<
" c_to_g handshake finished on coupler rank "
<<
worldRank
<<
std
::
endl
;
adios2
::
Engine
g_to_c_engine
=
g_to_c_io
.
Open
(
"g_to_c"
,
adios2
::
Mode
::
Read
);
std
::
cout
<<
" ======================== "
<<
std
::
endl
;
std
::
cout
<<
" g_to_c handshake finished on coupler rank "
<<
worldRank
<<
std
::
endl
;
for
(
int
i
=
0
;
i
<
steps
;
++
i
)
{
...
...
@@ -152,13 +137,7 @@ void xgc(const Dims &shape, const Dims &start, const Dims &count,
adios2
::
Engine
x_to_c_engine
=
x_to_c_io
.
Open
(
"x_to_c"
,
adios2
::
Mode
::
Write
);
std
::
cout
<<
" ======================== "
<<
std
::
endl
;
std
::
cout
<<
" x_to_c handshake finished on xgc rank "
<<
worldRank
<<
std
::
endl
;
adios2
::
Engine
c_to_x_engine
=
c_to_x_io
.
Open
(
"c_to_x"
,
adios2
::
Mode
::
Read
);
std
::
cout
<<
" ======================== "
<<
std
::
endl
;
std
::
cout
<<
" c_to_x handshake finished on xgc rank "
<<
worldRank
<<
std
::
endl
;
for
(
int
i
=
0
;
i
<
steps
;
++
i
)
{
...
...
@@ -168,13 +147,13 @@ void xgc(const Dims &shape, const Dims &start, const Dims &count,
x_to_c_engine
.
EndStep
();
c_to_x_engine
.
BeginStep
();
auto
c_to_x_var
=
x
_to_
c
_io
.
InquireVariable
<
float
>
(
"c_to_x"
);
c_to_x_data
.
resize
(
std
::
accumulate
(
c_to_x_var
.
Shape
()
.
begin
(),
c_to_x_var
.
Shape
()
.
end
(),
1
,
std
::
multiplies
<
size_t
>
()));
auto
c_to_x_var
=
c
_to_
x
_io
.
InquireVariable
<
float
>
(
"c_to_x"
);
auto
readShape
=
c_to_x_var
.
Shape
()
;
c_to_x_data
.
resize
(
std
::
accumulate
(
readShape
.
begin
(),
readShape
.
end
(),
1
,
std
::
multiplies
<
size_t
>
()));
c_to_x_engine
.
Get
(
c_to_x_var
,
c_to_x_data
.
data
(),
adios2
::
Mode
::
Sync
);
VerifyData
(
c_to_x_data
.
data
(),
i
,
Dims
(
c_to_x_var
.
Shape
()
.
size
(),
0
),
c_to_x_var
.
Shape
(),
c_to_x_var
.
Shape
()
,
mpiRank
);
VerifyData
(
c_to_x_data
.
data
(),
i
,
Dims
(
read
Shape
.
size
(),
0
),
readShape
,
read
Shape
,
mpiRank
);
c_to_x_engine
.
EndStep
();
}
...
...
@@ -200,14 +179,8 @@ void gene(const Dims &shape, const Dims &start, const Dims &count,
g_to_c_io
.
DefineVariable
<
float
>
(
"g_to_c"
,
shape
,
start
,
count
);
adios2
::
Engine
c_to_g_engine
=
c_to_g_io
.
Open
(
"c_to_g"
,
adios2
::
Mode
::
Read
);
std
::
cout
<<
" ======================== "
<<
std
::
endl
;
std
::
cout
<<
" c_to_g handshake finished on gene rank "
<<
worldRank
<<
std
::
endl
;
adios2
::
Engine
g_to_c_engine
=
g_to_c_io
.
Open
(
"g_to_c"
,
adios2
::
Mode
::
Write
);
std
::
cout
<<
" ======================== "
<<
std
::
endl
;
std
::
cout
<<
" g_to_c handshake finished on gene rank "
<<
worldRank
<<
std
::
endl
;
size_t
datasize
=
std
::
accumulate
(
shape
.
begin
(),
shape
.
end
(),
1
,
std
::
multiplies
<
size_t
>
());
...
...
@@ -218,12 +191,12 @@ void gene(const Dims &shape, const Dims &start, const Dims &count,
{
c_to_g_engine
.
BeginStep
(
StepMode
::
Read
,
5
);
auto
c_to_g_var
=
c_to_g_io
.
InquireVariable
<
float
>
(
"c_to_g"
);
c_to_g_data
.
resize
(
std
::
accumulate
(
c_to_g_var
.
Shape
()
.
begin
(),
c_to_g_var
.
Shape
()
.
end
(),
1
,
std
::
multiplies
<
size_t
>
()));
auto
readShape
=
c_to_g_var
.
Shape
()
;
c_to_g_data
.
resize
(
std
::
accumulate
(
readShape
.
begin
(),
readShape
.
end
(),
1
,
std
::
multiplies
<
size_t
>
()));
c_to_g_engine
.
Get
(
c_to_g_var
,
c_to_g_data
.
data
(),
adios2
::
Mode
::
Sync
);
VerifyData
(
c_to_g_data
.
data
(),
i
,
Dims
(
s
hape
.
size
(),
0
),
shape
,
s
hape
,
mpiRank
);
VerifyData
(
c_to_g_data
.
data
(),
i
,
Dims
(
readS
hape
.
size
(),
0
),
readS
hape
,
readShape
,
mpiRank
);
c_to_g_engine
.
EndStep
();
g_to_c_engine
.
BeginStep
();
...
...
@@ -277,8 +250,7 @@ TEST_F(SscEngineTest, TestSscXgc3Way)
adios2
::
Params
engineParams
=
{{
"RendezvousAppCount"
,
"2"
},
{
"MaxStreamsPerApp"
,
"4"
},
{
"OpenTimeoutSecs"
,
"3"
},
{
"Verbose"
,
"10"
}};
std
::
cout
<<
"Rank "
<<
worldRank
<<
" launched coupler"
<<
std
::
endl
;
{
"Verbose"
,
"0"
}};
coupler
(
shape
,
start
,
count
,
steps
,
engineParams
);
}
...
...
@@ -290,8 +262,7 @@ TEST_F(SscEngineTest, TestSscXgc3Way)
adios2
::
Params
engineParams
=
{{
"RendezvousAppCount"
,
"2"
},
{
"MaxStreamsPerApp"
,
"4"
},
{
"OpenTimeoutSecs"
,
"3"
},
{
"Verbose"
,
"10"
}};
std
::
cout
<<
"Rank "
<<
worldRank
<<
" launched gene"
<<
std
::
endl
;
{
"Verbose"
,
"0"
}};
gene
(
shape
,
start
,
shape
,
steps
,
engineParams
);
}
...
...
@@ -303,8 +274,7 @@ TEST_F(SscEngineTest, TestSscXgc3Way)
adios2
::
Params
engineParams
=
{{
"RendezvousAppCount"
,
"2"
},
{
"MaxStreamsPerApp"
,
"4"
},
{
"OpenTimeoutSecs"
,
"3"
},
{
"Verbose"
,
"10"
}};
std
::
cout
<<
"Rank "
<<
worldRank
<<
" launched xgc"
<<
std
::
endl
;
{
"Verbose"
,
"0"
}};
xgc
(
shape
,
start
,
count
,
steps
,
engineParams
);
}
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment