35#include <boost/algorithm/string/predicate.hpp>
36#include <boost/core/ignore_unused.hpp>
42#include <boost/algorithm/string/classification.hpp>
43#include <boost/algorithm/string/split.hpp>
55std::map<size_t, map<OpImpTimingKey, OperatorImpMap>>
61 : m_shapeDim(shapedim)
63 map<ElmtOrder, ImplementationType> defaults, defaultsPhysDeriv,
65 bool verbose = (pSession.get()) &&
66 (pSession->DefinesCmdLineArgument(
"verbose")) &&
67 (pSession->GetComm()->GetRank() == 0);
71 m_timeLevel = pSession.get() ? pSession->GetTimeLevel() : 0;
74 map<string, LibUtilities::ShapeType> elTypes;
84 for (
auto &it2 : elTypes)
93 for (
auto &it2 : elTypes)
103 map<string, OperatorType> opTypes;
120 map<string, ImplementationType> impTypes;
128 if ((defaultType ==
eNoImpType) && (pSession.get()))
130 TiXmlDocument &doc = pSession->GetDocument();
131 TiXmlHandle docHandle(&doc);
132 TiXmlElement *master = docHandle.FirstChildElement(
"NEKTAR").Element();
133 ASSERTL0(master,
"Unable to find NEKTAR tag in file.");
134 bool WriteFullCollections =
false;
136 TiXmlElement *xmlCol = master->FirstChildElement(
"COLLECTIONS");
144 const char *maxSize = xmlCol->Attribute(
"MAXSIZE");
147 const char *defaultImpl = xmlCol->Attribute(
"DEFAULT");
154 const std::string collinfo = string(defaultImpl);
155 m_autotune = boost::iequals(collinfo,
"auto");
159 bool collectionFound{
false};
168 collectionFound =
true;
174 "Unknown default collection scheme: " + collinfo);
178 for (
auto &it2 : elTypes)
189 const char *write = xmlCol->Attribute(
"WRITE");
190 if (write && boost::iequals(write,
"true"))
192 WriteFullCollections =
true;
201 if (WriteFullCollections)
208 for (
auto &eIt : mIt.second)
212 <<
" order " << eIt.first.second <<
" -> "
225 bool verboseHeader =
true;
226 map<string, LibUtilities::ShapeType> elTypes;
235 map<string, OperatorType> opTypes;
241 map<string, ImplementationType> impTypes;
247 TiXmlElement *elmt = xmlCol->FirstChildElement();
250 string tagname = elmt->ValueStr();
252 ASSERTL0(boost::iequals(tagname,
"OPERATOR"),
253 "Only OPERATOR tags are supported inside the "
256 const char *attr = elmt->Attribute(
"TYPE");
257 ASSERTL0(attr,
"Missing TYPE in OPERATOR tag.");
261 "Unknown OPERATOR type " + opType +
".");
265 TiXmlElement *elmt2 = elmt->FirstChildElement();
267 map<int, pair<int, std::string>> verboseWrite;
270 string tagname = elmt2->ValueStr();
271 ASSERTL0(boost::iequals(tagname,
"ELEMENT"),
272 "Only ELEMENT tags are supported inside the "
275 const char *attr = elmt2->Attribute(
"TYPE");
276 ASSERTL0(attr,
"Missing TYPE in ELEMENT tag.");
279 auto it2 = elTypes.find(elType);
280 ASSERTL0(it2 != elTypes.end(),
"Unknown element type " + elType +
284 const char *attr2 = elmt2->Attribute(
"IMPTYPE");
285 ASSERTL0(attr2,
"Missing IMPTYPE in ELEMENT tag.");
286 string impType(attr2);
287 ASSERTL0(impTypes.count(impType) > 0,
288 "Unknown IMPTYPE type " + impType +
".");
290 const char *attr3 = elmt2->Attribute(
"ORDER");
291 ASSERTL0(attr3,
"Missing ORDER in ELEMENT tag.");
299 global[ot][
ElmtOrder(it2->second, -1)] = impTypes[impType];
303 verboseWrite[it2->second] =
304 pair<int, std::string>(-1, impType);
309 vector<unsigned int> orders;
312 for (
int i = 0; i < orders.size(); ++i)
314 global[ot][
ElmtOrder(it2->second, orders[i])] =
319 verboseWrite[it2->second] =
320 pair<int, std::string>(orders[i], impType);
326 elmt2 = elmt2->NextSiblingElement();
329 if (verboseWrite.size())
333 cout <<
"Collection settings from file: " << endl;
334 verboseHeader =
false;
339 for (
auto &it : verboseWrite)
342 <<
" order " << it.second.first <<
" -> "
343 << it.second.second << endl;
347 elmt = elmt->NextSiblingElement();
355 ElmtOrder searchKey(pExp->DetShapeType(), pExp->GetBasisNumModes(0));
356 ElmtOrder defSearch(pExp->DetShapeType(), -1);
362 auto it2 = it.second.find(searchKey);
364 if (it2 == it.second.end())
366 it2 = it.second.find(defSearch);
367 if (it2 == it.second.end())
374 impType = it2->second;
379 impType = it2->second;
382 ret[it.first] = impType;
389 vector<StdRegions::StdExpansionSharedPtr> pCollExp,
392 boost::ignore_unused(impTypes);
408 pCollExp.size() * max(pExp->GetNcoeffs(), pExp->GetTotPoints());
418 cout <<
"Collection Implementation for "
420 for (
int i = 0; i < pExp->GetNumBases(); ++i)
422 cout << pExp->GetBasis(i)->GetNumModes() <<
" ";
425 <<
" for ngeoms = " << pCollExp.size() << endl;
440 OperatorKey opKey(pCollExp[0]->DetShapeType(), opType, impType,
441 pCollExp[0]->IsNodalNonTensorialExp());
445 impTypes[opType] = impType;
454 coll.push_back(collLoc);
465 coll[0].ApplyOperator(OpType, inarray, outarray1, outarray2, outarray3);
470 Ntest[i] = max((
int)(0.25 / oneTest), 1);
481 <<
"\t (IterLocExp, IterStdExp, "
482 "StdMat, SumFac, MatrixFree)"
492 for (
int imp = 0; imp < coll.size(); ++imp)
494 if (coll[imp].HasOperator(OpType))
497 for (
int n = 0; n < Ntest[i]; ++n)
499 coll[imp].ApplyOperator(OpType, inarray, outarray1,
500 outarray2, outarray3);
507 timing[imp] = 1000.0;
512 int minImp =
Vmath::Imin(coll.size(), timing, 1) + 1;
518 for (
int j = 0; j < coll.size(); ++j)
520 if (timing[j] > 999.0)
528 if (j != coll.size() - 1)
552 bool parallelInTime = comm->GetSize() != comm->GetSpaceComm()->GetSize();
554 if (parallelInTime && comm->GetTimeComm()->GetRank() > 0)
560 std::string outname = sessName.substr(0, sessName.find(
"_xml/")) +
".opt";
564 TiXmlElement *xmlCol =
new TiXmlElement(
"COLLECTIONS");
566 int rank = comm->GetSpaceComm()->GetRank();
569 if (!doc.LoadFile(outname))
571 TiXmlDeclaration *decl =
new TiXmlDeclaration(
"1.0",
"utf-8",
"");
572 doc.LinkEndChild(decl);
573 root =
new TiXmlElement(
"NEKTAR");
574 doc.LinkEndChild(root);
575 root->LinkEndChild(xmlCol);
579 xmlCol =
new TiXmlElement(
"TIMELEVEL");
580 xmlCol->SetAttribute(
"VALUE", 0);
581 root->FirstChildElement(
"COLLECTIONS")->LinkEndChild(xmlCol);
586 root = doc.FirstChildElement(
"NEKTAR");
587 xmlCol = root->FirstChildElement(
"COLLECTIONS");
593 bool verbose =
false;
599 xmlCol =
new TiXmlElement(
"TIMELEVEL");
601 root->FirstChildElement(
"COLLECTIONS")->LinkEndChild(xmlCol);
607 map<LibUtilities::ShapeType, int> ShapeMaxSize;
610 bool updateShape =
true;
615 if (ShapeMaxSize.count(shape))
617 int ngeoms = opimp.first.GetNGeoms();
618 if (ngeoms > ShapeMaxSize[shape])
620 ShapeMaxSize[shape] = ngeoms;
630 for (
auto &op : opimp.second)
632 global[op.first][
ElmtOrder(shape, -1)] = op.second;
638 for (
auto &op : global)
643 for (
auto &el : op.second)
645 ElmtImp[el.first.first] = el.second;
646 ElmtDef[el.first.first] =
true;
656 if ((ElmtImp[i] != -1) && (ElmtDef[i] ==
false))
671 map<LibUtilities::ShapeType, string> ShapeLetMap = {
680 for (
auto &op : global)
682 TiXmlElement *ColOp =
new TiXmlElement(
"OPERATOR");
683 xmlCol->LinkEndChild(ColOp);
686 for (
auto &el : op.second)
688 TiXmlElement *ElmtOp =
new TiXmlElement(
"ELEMENT");
689 ColOp->LinkEndChild(ElmtOp);
691 ElmtOp->SetAttribute(
"TYPE", ShapeLetMap[el.first.first]);
692 ElmtOp->SetAttribute(
"ORDER",
"*");
693 ElmtOp->SetAttribute(
"IMPTYPE",
698 doc.SaveFile(outname);
#define ASSERTL0(condition, msg)
COLLECTIONS_EXPORT void Initialise(const OperatorType opType, StdRegions::FactorMap factors=StdRegions::NullFactorMap)
COLLECTIONS_EXPORT OperatorImpMap SetWithTimings(std::vector< StdRegions::StdExpansionSharedPtr > pGeom, OperatorImpMap &impTypes, bool verbose=true)
ImplementationType m_defaultType
COLLECTIONS_EXPORT void UpdateOptFile(std::string sessName, LibUtilities::CommSharedPtr &comm)
static std::map< size_t, std::map< OpImpTimingKey, OperatorImpMap > > m_opImpMap
COLLECTIONS_EXPORT CollectionOptimisation(LibUtilities::SessionReaderSharedPtr pSession, const int shapedim, ImplementationType defaultType=eStdMat)
COLLECTIONS_EXPORT OperatorImpMap GetOperatorImpMap(StdRegions::StdExpansionSharedPtr pExp)
Get Operator Implementation Map from XMl or using default;.
std::map< OperatorType, std::map< ElmtOrder, ImplementationType > > GlobalOpMap
std::pair< LibUtilities::ShapeType, int > ElmtOrder
void ReadCollOps(TiXmlElement *xmlCol, GlobalOpMap &global, bool verbose)
static void GetXMLElementTimeLevel(TiXmlElement *&element, const size_t timeLevel, const bool disableCheck=true)
Get XML elment time level (Parallel-in-Time)
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
static bool GenerateSeqVector(const std::string &str, std::vector< unsigned int > &out)
Takes a comma-separated compressed string and converts it to entries in a vector.
std::map< OperatorType, ImplementationType > OperatorImpMap
const char *const ImplementationTypeMap[]
const char *const ImplementationTypeMap1[]
@ SIZE_ImplementationType
const char *const OperatorTypeMap[]
std::tuple< LibUtilities::ShapeType, OperatorType, ImplementationType, ExpansionIsNodal > OperatorKey
Key for describing an Operator.
OperatorFactory & GetOperatorFactory()
Returns the singleton Operator factory object.
const char *const OperatorTypeMap1[]
std::vector< Collection > CollectionVector
const char *const ShapeTypeMap[SIZE_ShapeType]
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
constexpr unsigned int ShapeTypeDimMap[SIZE_ShapeType]
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::map< ConstFactorType, NekDouble > ConstFactorMap
StdRegions::ConstFactorMap factors
The above copyright notice and this permission notice shall be included.
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.