35#include <boost/algorithm/string/predicate.hpp>
41#include <boost/algorithm/string/classification.hpp>
42#include <boost/algorithm/string/split.hpp>
52std::map<size_t, map<OpImpTimingKey, OperatorImpMap>>
58 : m_shapeDim(shapedim)
60 map<ElmtOrder, ImplementationType> defaults, defaultsPhysDeriv,
61 defaultsHelmholtz, defaultsPhysInterp1DScaled;
62 bool verbose = (pSession.get()) &&
63 (pSession->DefinesCmdLineArgument(
"verbose")) &&
64 (pSession->GetComm()->GetRank() == 0);
68 m_timeLevel = pSession.get() ? pSession->GetTimeLevel() : 0;
71 map<string, LibUtilities::ShapeType> elTypes;
81 for (
auto &it2 : elTypes)
91 for (
auto &it2 : elTypes)
100 defaultsPhysInterp1DScaled[
ElmtOrder(it2.second, -1)] =
105 map<string, OperatorType> opTypes;
125 map<string, ImplementationType> impTypes;
133 if ((defaultType ==
eNoImpType) && (pSession.get()))
135 TiXmlDocument &doc = pSession->GetDocument();
136 TiXmlHandle docHandle(&doc);
137 TiXmlElement *master = docHandle.FirstChildElement(
"NEKTAR").Element();
138 ASSERTL0(master,
"Unable to find NEKTAR tag in file.");
139 bool WriteFullCollections =
false;
141 TiXmlElement *xmlCol = master->FirstChildElement(
"COLLECTIONS");
149 const char *maxSize = xmlCol->Attribute(
"MAXSIZE");
152 const char *defaultImpl = xmlCol->Attribute(
"DEFAULT");
159 const std::string collinfo = string(defaultImpl);
160 m_autotune = boost::iequals(collinfo,
"auto");
164 bool collectionFound{
false};
173 collectionFound =
true;
179 "Unknown default collection scheme: " + collinfo);
183 for (
auto &it2 : elTypes)
194 const char *write = xmlCol->Attribute(
"WRITE");
195 if (write && boost::iequals(write,
"true"))
197 WriteFullCollections =
true;
206 if (WriteFullCollections)
213 for (
auto &eIt : mIt.second)
217 <<
" order " << eIt.first.second <<
" -> "
230 bool verboseHeader =
true;
231 map<string, LibUtilities::ShapeType> elTypes;
240 map<string, OperatorType> opTypes;
246 map<string, ImplementationType> impTypes;
252 TiXmlElement *elmt = xmlCol->FirstChildElement();
255 string tagname = elmt->ValueStr();
257 ASSERTL0(boost::iequals(tagname,
"OPERATOR"),
258 "Only OPERATOR tags are supported inside the "
261 const char *attr = elmt->Attribute(
"TYPE");
262 ASSERTL0(attr,
"Missing TYPE in OPERATOR tag.");
266 "Unknown OPERATOR type " + opType +
".");
270 TiXmlElement *elmt2 = elmt->FirstChildElement();
272 map<int, pair<int, std::string>> verboseWrite;
275 string tagname = elmt2->ValueStr();
276 ASSERTL0(boost::iequals(tagname,
"ELEMENT"),
277 "Only ELEMENT tags are supported inside the "
280 const char *attr = elmt2->Attribute(
"TYPE");
281 ASSERTL0(attr,
"Missing TYPE in ELEMENT tag.");
284 auto it2 = elTypes.find(elType);
285 ASSERTL0(it2 != elTypes.end(),
"Unknown element type " + elType +
289 const char *attr2 = elmt2->Attribute(
"IMPTYPE");
290 ASSERTL0(attr2,
"Missing IMPTYPE in ELEMENT tag.");
291 string impType(attr2);
292 ASSERTL0(impTypes.count(impType) > 0,
293 "Unknown IMPTYPE type " + impType +
".");
295 const char *attr3 = elmt2->Attribute(
"ORDER");
296 ASSERTL0(attr3,
"Missing ORDER in ELEMENT tag.");
304 global[ot][
ElmtOrder(it2->second, -1)] = impTypes[impType];
308 verboseWrite[it2->second] =
309 pair<int, std::string>(-1, impType);
314 vector<unsigned int> orders;
317 for (
int i = 0; i < orders.size(); ++i)
319 global[ot][
ElmtOrder(it2->second, orders[i])] =
324 verboseWrite[it2->second] =
325 pair<int, std::string>(orders[i], impType);
331 elmt2 = elmt2->NextSiblingElement();
334 if (verboseWrite.size())
338 cout <<
"Collection settings from file: " << endl;
339 verboseHeader =
false;
344 for (
auto &it : verboseWrite)
347 <<
" order " << it.second.first <<
" -> "
348 << it.second.second << endl;
352 elmt = elmt->NextSiblingElement();
360 ElmtOrder searchKey(pExp->DetShapeType(), pExp->GetBasisNumModes(0));
361 ElmtOrder defSearch(pExp->DetShapeType(), -1);
367 auto it2 = it.second.find(searchKey);
369 if (it2 == it.second.end())
371 it2 = it.second.find(defSearch);
372 if (it2 == it.second.end())
379 impType = it2->second;
384 impType = it2->second;
387 ret[it.first] = impType;
394 vector<StdRegions::StdExpansionSharedPtr> pCollExp,
414 cout <<
"Collection Implementation for "
416 for (
int i = 0; i < pExp->GetNumBases(); ++i)
418 cout << pExp->GetBasis(i)->GetNumModes() <<
" ";
421 <<
" for ngeoms = " << pCollExp.size() << endl;
432 int qpInsideElmt_idir;
434 int newQpInsideElmt_idir;
436 int newTotQpInsideElmt{1};
438 for (
int i = 0; i < pExp->GetShapeDimension(); ++i)
440 qpInsideElmt_idir = pExp->GetNumPoints(i);
441 newQpInsideElmt_idir =
443 newTotQpInsideElmt *= newQpInsideElmt_idir;
445 int maxsize = pCollExp.size() * max(pExp->GetNcoeffs(), newTotQpInsideElmt);
456 for (
int i = 0; i < pExp->GetShapeDimension(); i++)
465 std::vector<OperatorType> opTypes;
469 OperatorKey opKey(pCollExp[0]->DetShapeType(), opType, impType,
470 pCollExp[0]->IsNodalNonTensorialExp());
474 impTypes[opType] = impType;
475 opTypes.push_back(opType);
480 for (
auto opType : opTypes)
489 coll.push_back(collLoc);
500 coll[0].ApplyOperator(OpType, inarray, outarray1, outarray2, outarray3);
505 Ntest[i] = max((
int)(0.25 / oneTest), 1);
516 <<
"\t (IterLocExp, IterStdExp, "
517 "StdMat, SumFac, MatrixFree)"
527 for (
int imp = 0; imp < coll.size(); ++imp)
529 if (coll[imp].HasOperator(OpType))
532 for (
int n = 0; n < Ntest[i]; ++n)
534 coll[imp].ApplyOperator(OpType, inarray, outarray1,
535 outarray2, outarray3);
542 timing[imp] = 1000.0;
547 int minImp =
Vmath::Imin(coll.size(), timing, 1) + 1;
553 for (
int j = 0; j < coll.size(); ++j)
555 if (timing[j] > 999.0)
563 if (j != coll.size() - 1)
587 if (comm->IsParallelInTime() && comm->GetTimeComm()->GetRank() > 0)
593 std::string outname = sessName.substr(0, sessName.find(
"_xml/")) +
".opt";
597 TiXmlElement *xmlCol =
new TiXmlElement(
"COLLECTIONS");
599 int rank = comm->GetSpaceComm()->GetRank();
602 if (!doc.LoadFile(outname))
604 TiXmlDeclaration *decl =
new TiXmlDeclaration(
"1.0",
"utf-8",
"");
605 doc.LinkEndChild(decl);
606 root =
new TiXmlElement(
"NEKTAR");
607 doc.LinkEndChild(root);
608 root->LinkEndChild(xmlCol);
609 if (comm->IsParallelInTime())
612 xmlCol =
new TiXmlElement(
"TIMELEVEL");
613 xmlCol->SetAttribute(
"VALUE", 0);
614 root->FirstChildElement(
"COLLECTIONS")->LinkEndChild(xmlCol);
619 root = doc.FirstChildElement(
"NEKTAR");
620 xmlCol = root->FirstChildElement(
"COLLECTIONS");
626 bool verbose =
false;
632 xmlCol =
new TiXmlElement(
"TIMELEVEL");
634 root->FirstChildElement(
"COLLECTIONS")->LinkEndChild(xmlCol);
640 map<LibUtilities::ShapeType, int> ShapeMaxSize;
643 bool updateShape =
true;
648 if (ShapeMaxSize.count(shape))
650 int ngeoms = opimp.first.GetNGeoms();
651 if (ngeoms > ShapeMaxSize[shape])
653 ShapeMaxSize[shape] = ngeoms;
663 for (
auto &op : opimp.second)
665 global[op.first][
ElmtOrder(shape, -1)] = op.second;
671 for (
auto &op : global)
676 for (
auto &el : op.second)
678 ElmtImp[el.first.first] = el.second;
679 ElmtDef[el.first.first] =
true;
689 if ((ElmtImp[i] != -1) && (ElmtDef[i] ==
false))
704 map<LibUtilities::ShapeType, string> ShapeLetMap = {
713 for (
auto &op : global)
715 TiXmlElement *ColOp =
new TiXmlElement(
"OPERATOR");
716 xmlCol->LinkEndChild(ColOp);
719 for (
auto &el : op.second)
721 TiXmlElement *ElmtOp =
new TiXmlElement(
"ELEMENT");
722 ColOp->LinkEndChild(ElmtOp);
724 ElmtOp->SetAttribute(
"TYPE", ShapeLetMap[el.first.first]);
725 ElmtOp->SetAttribute(
"ORDER",
"*");
726 ElmtOp->SetAttribute(
"IMPTYPE",
731 doc.SaveFile(outname);
#define ASSERTL0(condition, msg)
COLLECTIONS_EXPORT void UpdateVarcoeffs(const OperatorType opType, StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
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 enableCheck=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
@ eLinearAdvectionDiffusionReaction
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
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
StdRegions::ConstFactorMap factors
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.