Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Query fix #2758

Merged
merged 3 commits into from
Jun 9, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
83 changes: 77 additions & 6 deletions source/adios2/toolkit/query/Query.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -38,6 +38,9 @@ void QueryBase::ApplyOutputRegion(std::vector<Box<Dims>> &touchedBlocks,
}
bool QueryComposite::AddNode(QueryBase *var)
{
if (nullptr == var)
return false;

if (adios2::query::Relation::NOT == m_Relation)
{
// if (m_Nodes.size() > 0) return false;
Expand All @@ -54,6 +57,68 @@ void QueryComposite::BlockIndexEvaluate(adios2::core::IO &io,
adios2::core::Engine &reader,
std::vector<Box<Dims>> &touchedBlocks)
{
auto lf_ApplyAND = [&](std::vector<Box<Dims>> &touched,
const std::vector<Box<Dims>> &curr) -> void {
if (curr.size() == 0)
{
touched.clear();
return;
}

for (auto i = touched.size(); i >= 1; i--)
{
bool intersects = false;
for (auto b : curr)
{
adios2::Box<Dims> curr = GetIntersection(touched[i - 1], b);
if (curr.first.size() != 0) // has intersection
{
intersects = true;
break;
}
}
if (!intersects)
// it = touched.erase(it);
touched.erase(touched.begin() + i - 1);
// if (touched.end() == it)
// break;
}

for (auto b : curr)
{
for (auto it = touched.begin(); it != touched.end(); it++)
{
adios2::Box<Dims> curr = GetIntersection(*it, b);
if (curr.first.size() != 0) // has intersection
{
*it = curr;
}
}
}
}; // lf_ApplyAND

auto lf_ApplyOR = [&](std::vector<Box<Dims>> &touched,
const std::vector<Box<Dims>> &curr) -> void {
if (curr.size() == 0)
return;

for (auto b : curr)
{
bool duplicated = false;
for (auto box : touched)
{
if (adios2::helper::IdenticalBoxes(box, b))
{
duplicated = true;
continue;
}
}
if (!duplicated)
touched.push_back(b);
}
}; // lf_ApplyOR

/*
auto lf_ApplyRelation = [&](std::vector<Box<Dims>> &collection,
const Box<Dims> &block) -> void {
if (adios2::query::Relation::AND == m_Relation)
Expand All @@ -64,10 +129,15 @@ void QueryComposite::BlockIndexEvaluate(adios2::core::IO &io,
adios2::Box<Dims> curr = GetIntersection(*it, block);
// adios2::helper::IntersectionBox(*it, block);
if (curr.first.size() == 0) // no intersection
touchedBlocks.erase(it);
{
it = touchedBlocks.erase(it);
if (touchedBlocks.end() == it)
return;
}
else
*it = curr;
*it = curr;
}

return;
}

Expand All @@ -82,6 +152,7 @@ void QueryComposite::BlockIndexEvaluate(adios2::core::IO &io,
return;
}
}; // local
*/

if (m_Nodes.size() == 0)
return;
Expand Down Expand Up @@ -109,10 +180,10 @@ void QueryComposite::BlockIndexEvaluate(adios2::core::IO &io,
continue;
}

for (auto block : currBlocks)
{
lf_ApplyRelation(touchedBlocks, block);
}
if (adios2::query::Relation::AND == m_Relation)
lf_ApplyAND(touchedBlocks, currBlocks);
else if (adios2::query::Relation::OR == m_Relation)
lf_ApplyOR(touchedBlocks, currBlocks);
}
// plan to shift all var results to regions start at 0, and find out the
// overlapped regions boxes can be different size especially if they are
Expand Down
2 changes: 1 addition & 1 deletion source/adios2/toolkit/query/Query.h
Original file line number Diff line number Diff line change
Expand Up @@ -270,7 +270,7 @@ class QueryComposite : public QueryBase
bool IsCompatible(const adios2::Box<adios2::Dims> &box)
{
if (m_Nodes.size() == 0)
return false;
return true;
return (m_Nodes[0])->IsCompatible(box);
}

Expand Down
5 changes: 4 additions & 1 deletion source/adios2/toolkit/query/XmlWorker.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -94,6 +94,9 @@ void XmlWorker::ParseIONode(const pugi::xml_node &ioNode)
const pugi::xml_node &variable = qTagNode.child("var");
QueryVar *q =
ParseVarNode(variable, m_SourceReader->m_IO, *m_SourceReader);
if (!q)
continue;

if (ref.first.size() == 0)
{
ref = q->m_Selection;
Expand Down Expand Up @@ -140,7 +143,7 @@ QueryVar *XmlWorker::ParseVarNode(const pugi::xml_node &node,
if (varType == DataType::None)
{
std::cerr << "No such variable: " << variableName << std::endl;
return nullptr;
throw std::ios_base::failure("No such variable: " + variableName);
}
#define declare_type(T) \
if (varType == helper::GetDataType<T>()) \
Expand Down