Nektar++
Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
Nektar::Collections::CollectionOptimisation Class Reference

#include <CollectionOptimisation.h>

Public Member Functions

COLLECTIONS_EXPORT CollectionOptimisation (LibUtilities::SessionReaderSharedPtr pSession, const int shapedim, ImplementationType defaultType=eStdMat)
 
 ~CollectionOptimisation ()
 
ImplementationType GetDefaultImplementationType ()
 
size_t GetMaxCollectionSize ()
 
bool IsUsingAutotuning ()
 
COLLECTIONS_EXPORT OperatorImpMap GetOperatorImpMap (StdRegions::StdExpansionSharedPtr pExp)
 Get Operator Implementation Map from XMl or using default;. More...
 
COLLECTIONS_EXPORT OperatorImpMap SetWithTimings (std::vector< StdRegions::StdExpansionSharedPtr > pGeom, OperatorImpMap &impTypes, bool verbose=true)
 
COLLECTIONS_EXPORT void UpdateOptFile (std::string sessName, LibUtilities::CommSharedPtr &comm)
 

Private Types

typedef std::pair< LibUtilities::ShapeType, int > ElmtOrder
 
typedef std::map< OperatorType, std::map< ElmtOrder, ImplementationType > > GlobalOpMap
 

Private Member Functions

void ReadCollOps (TiXmlElement *xmlCol, GlobalOpMap &global, bool verbose)
 

Private Attributes

GlobalOpMap m_global
 
bool m_autotune
 
ImplementationType m_defaultType
 
size_t m_maxCollSize
 
size_t m_shapeDim
 
size_t m_timeLevel
 

Static Private Attributes

static std::map< size_t, std::map< OpImpTimingKey, OperatorImpMap > > m_opImpMap
 

Detailed Description

Definition at line 123 of file CollectionOptimisation.h.

Member Typedef Documentation

◆ ElmtOrder

Definition at line 162 of file CollectionOptimisation.h.

◆ GlobalOpMap

Definition at line 164 of file CollectionOptimisation.h.

Constructor & Destructor Documentation

◆ CollectionOptimisation()

Nektar::Collections::CollectionOptimisation::CollectionOptimisation ( LibUtilities::SessionReaderSharedPtr  pSession,
const int  shapedim,
ImplementationType  defaultType = eStdMat 
)

Definition at line 58 of file CollectionOptimisation.cpp.

61 : m_shapeDim(shapedim)
62{
63 map<ElmtOrder, ImplementationType> defaults, defaultsPhysDeriv,
64 defaultsHelmholtz;
65 bool verbose = (pSession.get()) &&
66 (pSession->DefinesCmdLineArgument("verbose")) &&
67 (pSession->GetComm()->GetRank() == 0);
68
69 m_autotune = false;
70 m_maxCollSize = 0;
71 m_timeLevel = pSession.get() ? pSession->GetTimeLevel() : 0;
72 m_defaultType = defaultType == eNoImpType ? eIterPerExp : defaultType;
73
74 map<string, LibUtilities::ShapeType> elTypes;
75 elTypes["S"] = LibUtilities::eSegment;
76 elTypes["T"] = LibUtilities::eTriangle;
77 elTypes["Q"] = LibUtilities::eQuadrilateral;
78 elTypes["A"] = LibUtilities::eTetrahedron;
79 elTypes["P"] = LibUtilities::ePyramid;
80 elTypes["R"] = LibUtilities::ePrism;
81 elTypes["H"] = LibUtilities::eHexahedron;
82
83 // Set defaults for all element types.
84 for (auto &it2 : elTypes)
85 {
86 defaults[ElmtOrder(it2.second, -1)] = m_defaultType;
87 defaultsPhysDeriv[ElmtOrder(it2.second, -1)] = m_defaultType;
88 defaultsHelmholtz[ElmtOrder(it2.second, -1)] = m_defaultType;
89 }
90
91 if (defaultType == eNoImpType)
92 {
93 for (auto &it2 : elTypes)
94 {
95 // use Nocollection for Phys Deriv
96 defaultsPhysDeriv[ElmtOrder(it2.second, -1)] = eNoCollection;
97
98 // Use IterPerExp
99 defaultsHelmholtz[ElmtOrder(it2.second, -1)] = eMatrixFree;
100 }
101 }
102
103 map<string, OperatorType> opTypes;
104 for (int i = 0; i < SIZE_OperatorType; ++i)
105 {
106 opTypes[OperatorTypeMap[i]] = (OperatorType)i;
107 switch ((OperatorType)i)
108 {
109 case eHelmholtz:
110 m_global[(OperatorType)i] = defaultsHelmholtz;
111 break;
112 case ePhysDeriv:
113 m_global[(OperatorType)i] = defaultsPhysDeriv;
114 break;
115 default:
116 m_global[(OperatorType)i] = defaults;
117 }
118 }
119
120 map<string, ImplementationType> impTypes;
121 for (int i = 0; i < SIZE_ImplementationType; ++i)
122 {
124 }
125
126 // turn off file reader if dummy pointer is given or if default
127 // option is passed and by default calling argument.
128 if ((defaultType == eNoImpType) && (pSession.get()))
129 {
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;
135
136 TiXmlElement *xmlCol = master->FirstChildElement("COLLECTIONS");
139
140 // Check if user has specified some options
141 if (xmlCol)
142 {
143 // Set the maxsize and default implementation type if provided
144 const char *maxSize = xmlCol->Attribute("MAXSIZE");
145 m_maxCollSize = (maxSize ? atoi(maxSize) : 0);
146
147 const char *defaultImpl = xmlCol->Attribute("DEFAULT");
148 m_defaultType = defaultType;
149
150 // If user has specified a default impl type or autotuning
151 // and set this default across all operators.
152 if (defaultImpl)
153 {
154 const std::string collinfo = string(defaultImpl);
155 m_autotune = boost::iequals(collinfo, "auto");
156
157 if (!m_autotune)
158 {
159 bool collectionFound{false};
160 for (int i = 1; i < Collections::SIZE_ImplementationType;
161 ++i)
162 {
163 if (boost::iequals(
164 collinfo,
166 {
168 collectionFound = true;
169 break;
170 }
171 }
172
173 ASSERTL0(collectionFound,
174 "Unknown default collection scheme: " + collinfo);
175
176 defaults.clear();
177 // Override default types
178 for (auto &it2 : elTypes)
179 {
180 defaults[ElmtOrder(it2.second, -1)] = m_defaultType;
181 }
182
183 for (int i = 0; i < SIZE_OperatorType; ++i)
184 {
185 m_global[(OperatorType)i] = defaults;
186 }
187 }
188 }
189 const char *write = xmlCol->Attribute("WRITE");
190 if (write && boost::iequals(write, "true"))
191 {
192 WriteFullCollections = true;
193 }
194
195 // Now process operator-specific implementation selections
196 ReadCollOps(xmlCol, m_global, verbose);
197
198 // Print out operator map
199 if (verbose)
200 {
201 if (WriteFullCollections)
202 {
203 for (auto &mIt : m_global)
204 {
205 cout << "Operator " << OperatorTypeMap[mIt.first] << ":"
206 << endl;
207
208 for (auto &eIt : mIt.second)
209 {
210 cout << "- "
211 << LibUtilities::ShapeTypeMap[eIt.first.first]
212 << " order " << eIt.first.second << " -> "
213 << ImplementationTypeMap[eIt.second] << endl;
214 }
215 }
216 }
217 }
218 }
219 }
220}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
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)
const char *const ImplementationTypeMap[]
Definition: Operator.h:94
const char *const OperatorTypeMap[]
Definition: Operator.h:76
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:79

References ASSERTL0, Nektar::Collections::eHelmholtz, Nektar::LibUtilities::eHexahedron, Nektar::Collections::eIterPerExp, Nektar::Collections::eMatrixFree, Nektar::Collections::eNoCollection, Nektar::Collections::eNoImpType, Nektar::Collections::ePhysDeriv, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::ePyramid, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eSegment, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, Nektar::LibUtilities::SessionReader::GetXMLElementTimeLevel(), Nektar::Collections::ImplementationTypeMap, m_autotune, m_defaultType, m_global, m_maxCollSize, m_timeLevel, Nektar::Collections::OperatorTypeMap, ReadCollOps(), Nektar::LibUtilities::ShapeTypeMap, Nektar::Collections::SIZE_ImplementationType, and Nektar::Collections::SIZE_OperatorType.

◆ ~CollectionOptimisation()

Nektar::Collections::CollectionOptimisation::~CollectionOptimisation ( )
inline

Definition at line 131 of file CollectionOptimisation.h.

131{};

Member Function Documentation

◆ GetDefaultImplementationType()

ImplementationType Nektar::Collections::CollectionOptimisation::GetDefaultImplementationType ( )
inline

Definition at line 133 of file CollectionOptimisation.h.

134 {
135 return m_defaultType;
136 }

References m_defaultType.

◆ GetMaxCollectionSize()

size_t Nektar::Collections::CollectionOptimisation::GetMaxCollectionSize ( )
inline

Definition at line 138 of file CollectionOptimisation.h.

139 {
140 return m_maxCollSize;
141 }

References m_maxCollSize.

Referenced by Nektar::MultiRegions::ExpList::CreateCollections().

◆ GetOperatorImpMap()

OperatorImpMap Nektar::Collections::CollectionOptimisation::GetOperatorImpMap ( StdRegions::StdExpansionSharedPtr  pExp)

Get Operator Implementation Map from XMl or using default;.

Definition at line 351 of file CollectionOptimisation.cpp.

353{
354 OperatorImpMap ret;
355 ElmtOrder searchKey(pExp->DetShapeType(), pExp->GetBasisNumModes(0));
356 ElmtOrder defSearch(pExp->DetShapeType(), -1);
357
358 for (auto &it : m_global)
359 {
360 ImplementationType impType;
361
362 auto it2 = it.second.find(searchKey);
363
364 if (it2 == it.second.end())
365 {
366 it2 = it.second.find(defSearch);
367 if (it2 == it.second.end())
368 {
369 // Shouldn't be able to reach here.
370 impType = eNoCollection;
371 }
372 else
373 {
374 impType = it2->second;
375 }
376 }
377 else
378 {
379 impType = it2->second;
380 }
381
382 ret[it.first] = impType;
383 }
384
385 return ret;
386}
std::map< OperatorType, ImplementationType > OperatorImpMap
Definition: Operator.h:112

References Nektar::Collections::eNoCollection, and m_global.

Referenced by Nektar::HexCollectionTests::BOOST_AUTO_TEST_CASE(), Nektar::PrismCollectionTests::BOOST_AUTO_TEST_CASE(), Nektar::PyrCollectionTests::BOOST_AUTO_TEST_CASE(), Nektar::QuadCollectionTests::BOOST_AUTO_TEST_CASE(), Nektar::SegCollectionTests::BOOST_AUTO_TEST_CASE(), Nektar::TetCollectionTests::BOOST_AUTO_TEST_CASE(), Nektar::TriCollectionTests::BOOST_AUTO_TEST_CASE(), and Nektar::MultiRegions::ExpList::CreateCollections().

◆ IsUsingAutotuning()

bool Nektar::Collections::CollectionOptimisation::IsUsingAutotuning ( )
inline

Definition at line 143 of file CollectionOptimisation.h.

144 {
145 return m_autotune;
146 }

References m_autotune.

Referenced by Nektar::MultiRegions::ExpList::CreateCollections().

◆ ReadCollOps()

void Nektar::Collections::CollectionOptimisation::ReadCollOps ( TiXmlElement *  xmlCol,
GlobalOpMap global,
bool  verbose 
)
private

Definition at line 222 of file CollectionOptimisation.cpp.

224{
225 bool verboseHeader = true;
226 map<string, LibUtilities::ShapeType> elTypes;
227 elTypes["S"] = LibUtilities::eSegment;
228 elTypes["T"] = LibUtilities::eTriangle;
229 elTypes["Q"] = LibUtilities::eQuadrilateral;
230 elTypes["A"] = LibUtilities::eTetrahedron;
231 elTypes["P"] = LibUtilities::ePyramid;
232 elTypes["R"] = LibUtilities::ePrism;
233 elTypes["H"] = LibUtilities::eHexahedron;
234
235 map<string, OperatorType> opTypes;
236 for (int i = 0; i < SIZE_OperatorType; ++i)
237 {
238 opTypes[OperatorTypeMap[i]] = (OperatorType)i;
239 }
240
241 map<string, ImplementationType> impTypes;
242 for (int i = 0; i < SIZE_ImplementationType; ++i)
243 {
245 }
246
247 TiXmlElement *elmt = xmlCol->FirstChildElement();
248 while (elmt)
249 {
250 string tagname = elmt->ValueStr();
251
252 ASSERTL0(boost::iequals(tagname, "OPERATOR"),
253 "Only OPERATOR tags are supported inside the "
254 "COLLECTIONS tag.");
255
256 const char *attr = elmt->Attribute("TYPE");
257 ASSERTL0(attr, "Missing TYPE in OPERATOR tag.");
258 string opType(attr);
259
260 ASSERTL0(opTypes.count(opType) > 0,
261 "Unknown OPERATOR type " + opType + ".");
262
263 OperatorType ot = opTypes[opType];
264
265 TiXmlElement *elmt2 = elmt->FirstChildElement();
266
267 map<int, pair<int, std::string>> verboseWrite;
268 while (elmt2)
269 {
270 string tagname = elmt2->ValueStr();
271 ASSERTL0(boost::iequals(tagname, "ELEMENT"),
272 "Only ELEMENT tags are supported inside the "
273 "OPERATOR tag.");
274
275 const char *attr = elmt2->Attribute("TYPE");
276 ASSERTL0(attr, "Missing TYPE in ELEMENT tag.");
277
278 string elType(attr);
279 auto it2 = elTypes.find(elType);
280 ASSERTL0(it2 != elTypes.end(), "Unknown element type " + elType +
281 " in ELEMENT "
282 "tag");
283
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 + ".");
289
290 const char *attr3 = elmt2->Attribute("ORDER");
291 ASSERTL0(attr3, "Missing ORDER in ELEMENT tag.");
292 string order(attr3);
293
294 // load details relevant to this shape dimension.
295 if (LibUtilities::ShapeTypeDimMap[it2->second] == m_shapeDim)
296 {
297 if (order == "*")
298 {
299 global[ot][ElmtOrder(it2->second, -1)] = impTypes[impType];
300
301 if (verbose)
302 {
303 verboseWrite[it2->second] =
304 pair<int, std::string>(-1, impType);
305 }
306 }
307 else
308 {
309 vector<unsigned int> orders;
310 ParseUtils::GenerateSeqVector(order, orders);
311
312 for (int i = 0; i < orders.size(); ++i)
313 {
314 global[ot][ElmtOrder(it2->second, orders[i])] =
315 impTypes[impType];
316
317 if (verbose)
318 {
319 verboseWrite[it2->second] =
320 pair<int, std::string>(orders[i], impType);
321 }
322 }
323 }
324 }
325
326 elmt2 = elmt2->NextSiblingElement();
327 }
328
329 if (verboseWrite.size())
330 {
331 if (verboseHeader)
332 {
333 cout << "Collection settings from file: " << endl;
334 verboseHeader = false;
335 }
336
337 cout << "\t Operator " << OperatorTypeMap[ot] << ":" << endl;
338
339 for (auto &it : verboseWrite)
340 {
341 cout << "\t - " << LibUtilities::ShapeTypeMap[it.first]
342 << " order " << it.second.first << " -> "
343 << it.second.second << endl;
344 }
345 }
346
347 elmt = elmt->NextSiblingElement();
348 }
349}
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.
Definition: ParseUtils.cpp:105
constexpr unsigned int ShapeTypeDimMap[SIZE_ShapeType]
Definition: ShapeType.hpp:85

References ASSERTL0, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::ePyramid, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eSegment, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, Nektar::ParseUtils::GenerateSeqVector(), Nektar::Collections::ImplementationTypeMap, m_shapeDim, Nektar::Collections::OperatorTypeMap, Nektar::LibUtilities::ShapeTypeDimMap, Nektar::LibUtilities::ShapeTypeMap, Nektar::Collections::SIZE_ImplementationType, and Nektar::Collections::SIZE_OperatorType.

Referenced by CollectionOptimisation(), and UpdateOptFile().

◆ SetWithTimings()

OperatorImpMap Nektar::Collections::CollectionOptimisation::SetWithTimings ( std::vector< StdRegions::StdExpansionSharedPtr pGeom,
OperatorImpMap impTypes,
bool  verbose = true 
)

Definition at line 388 of file CollectionOptimisation.cpp.

391{
392 boost::ignore_unused(impTypes);
393
394 OperatorImpMap ret;
395
396 StdRegions::StdExpansionSharedPtr pExp = pCollExp[0];
397
398 // check to see if already defined for this expansion
399 OpImpTimingKey OpKey(pExp, pCollExp.size(), pExp->GetNumBases());
400 if (m_opImpMap.count(m_timeLevel) != 0 &&
401 m_opImpMap[m_timeLevel].count(OpKey) != 0)
402 {
403 ret = m_opImpMap[m_timeLevel][OpKey];
404 return ret;
405 }
406
407 int maxsize =
408 pCollExp.size() * max(pExp->GetNcoeffs(), pExp->GetTotPoints());
409 Array<OneD, NekDouble> inarray(maxsize, 1.0);
410 Array<OneD, NekDouble> outarray1(maxsize);
411 Array<OneD, NekDouble> outarray2(maxsize);
412 Array<OneD, NekDouble> outarray3(maxsize);
413
414 LibUtilities::Timer t;
415
416 if (verbose)
417 {
418 cout << "Collection Implementation for "
419 << LibUtilities::ShapeTypeMap[pExp->DetShapeType()] << " ( ";
420 for (int i = 0; i < pExp->GetNumBases(); ++i)
421 {
422 cout << pExp->GetBasis(i)->GetNumModes() << " ";
423 }
424 cout << ")"
425 << " for ngeoms = " << pCollExp.size() << endl;
426 }
427 // set up an array of collections
428 CollectionVector coll;
429
430 StdRegions::ConstFactorMap factors; // required for helmholtz operator
432
433 for (int imp = 1; imp < SIZE_ImplementationType; ++imp)
434 {
436 OperatorImpMap impTypes;
437 for (int i = 0; i < SIZE_OperatorType; ++i)
438 {
439 OperatorType opType = (OperatorType)i;
440 OperatorKey opKey(pCollExp[0]->DetShapeType(), opType, impType,
441 pCollExp[0]->IsNodalNonTensorialExp());
442
443 if (GetOperatorFactory().ModuleExists(opKey))
444 {
445 impTypes[opType] = impType;
446 }
447 }
448
449 Collection collLoc(pCollExp, impTypes);
450 for (int i = 0; i < SIZE_OperatorType; ++i)
451 {
452 collLoc.Initialise((OperatorType)i, factors);
453 }
454 coll.push_back(collLoc);
455 }
456
457 // Determine the number of tests to do in one second
458 Array<OneD, int> Ntest(SIZE_OperatorType);
459 for (int i = 0; i < SIZE_OperatorType; ++i)
460 {
461 OperatorType OpType = (OperatorType)i;
462
463 t.Start();
464
465 coll[0].ApplyOperator(OpType, inarray, outarray1, outarray2, outarray3);
466 t.Stop();
467
468 NekDouble oneTest = t.TimePerTest(1);
469
470 Ntest[i] = max((int)(0.25 / oneTest), 1);
471 }
472
473 Array<OneD, NekDouble> timing(SIZE_ImplementationType);
474
475 if (verbose)
476 {
477 cout << "\t "
478 << " Op. "
479 << ":\t"
480 << "opt. Impl."
481 << "\t (IterLocExp, IterStdExp, "
482 "StdMat, SumFac, MatrixFree)"
483 << endl;
484 }
485
486 // loop over all operators and determine fastest implementation
487 for (int i = 0; i < SIZE_OperatorType; ++i)
488 {
489 OperatorType OpType = (OperatorType)i;
490
491 // call collection implementation in thorugh ExpList.
492 for (int imp = 0; imp < coll.size(); ++imp)
493 {
494 if (coll[imp].HasOperator(OpType))
495 {
496 t.Start();
497 for (int n = 0; n < Ntest[i]; ++n)
498 {
499 coll[imp].ApplyOperator(OpType, inarray, outarray1,
500 outarray2, outarray3);
501 }
502 t.Stop();
503 timing[imp] = t.TimePerTest(Ntest[i]);
504 }
505 else
506 {
507 timing[imp] = 1000.0;
508 }
509 }
510 // determine optimal implementation. Note +1 to
511 // remove NoImplementationType flag
512 int minImp = Vmath::Imin(coll.size(), timing, 1) + 1;
513
514 if (verbose)
515 {
516 cout << "\t " << OperatorTypeMap1[i] << ": \t"
517 << ImplementationTypeMap1[minImp] << "\t (";
518 for (int j = 0; j < coll.size(); ++j)
519 {
520 if (timing[j] > 999.0)
521 {
522 cout << " -- ";
523 }
524 else
525 {
526 cout << timing[j];
527 }
528 if (j != coll.size() - 1)
529 {
530 cout << ", ";
531 }
532 }
533 cout << ")" << endl;
534 }
535
536 // set up new map
537 ret[OpType] = (ImplementationType)minImp;
538 }
539
540 // store map for use by another expansion.
541 if (m_opImpMap.count(m_timeLevel) == 0)
542 {
543 m_opImpMap[m_timeLevel] = map<OpImpTimingKey, OperatorImpMap>();
544 }
545 m_opImpMap[m_timeLevel][OpKey] = ret;
546 return ret;
547}
static std::map< size_t, std::map< OpImpTimingKey, OperatorImpMap > > m_opImpMap
const char *const ImplementationTypeMap1[]
Definition: Operator.h:101
std::tuple< LibUtilities::ShapeType, OperatorType, ImplementationType, ExpansionIsNodal > OperatorKey
Key for describing an Operator.
Definition: Operator.h:177
OperatorFactory & GetOperatorFactory()
Returns the singleton Operator factory object.
Definition: Operator.cpp:117
const char *const OperatorTypeMap1[]
Definition: Operator.h:80
std::vector< Collection > CollectionVector
Definition: Collection.h:110
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
StdRegions::ConstFactorMap factors
double NekDouble
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
Definition: Vmath.cpp:1018

References Nektar::StdRegions::eFactorLambda, Nektar::VarcoeffHashingTest::factors, Nektar::Collections::GetOperatorFactory(), Vmath::Imin(), Nektar::Collections::ImplementationTypeMap1, Nektar::Collections::Collection::Initialise(), m_opImpMap, m_timeLevel, Nektar::Collections::OperatorTypeMap1, Nektar::LibUtilities::ShapeTypeMap, Nektar::Collections::SIZE_ImplementationType, Nektar::Collections::SIZE_OperatorType, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Nektar::LibUtilities::Timer::TimePerTest().

Referenced by Nektar::MultiRegions::ExpList::CreateCollections().

◆ UpdateOptFile()

void Nektar::Collections::CollectionOptimisation::UpdateOptFile ( std::string  sessName,
LibUtilities::CommSharedPtr comm 
)

Definition at line 549 of file CollectionOptimisation.cpp.

551{
552 bool parallelInTime = comm->GetSize() != comm->GetSpaceComm()->GetSize();
553
554 if (parallelInTime && comm->GetTimeComm()->GetRank() > 0)
555 {
556 // No need to repeatly update the optfile for each time chunk.
557 return;
558 }
559
560 std::string outname = sessName.substr(0, sessName.find("_xml/")) + ".opt";
561
562 TiXmlDocument doc;
563 TiXmlElement *root;
564 TiXmlElement *xmlCol = new TiXmlElement("COLLECTIONS");
565 GlobalOpMap global;
566 int rank = comm->GetSpaceComm()->GetRank();
567 if (rank == 0)
568 {
569 if (!doc.LoadFile(outname)) // set up new file
570 {
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);
576 if (parallelInTime)
577 {
578 // Add timelevel tag
579 xmlCol = new TiXmlElement("TIMELEVEL");
580 xmlCol->SetAttribute("VALUE", 0);
581 root->FirstChildElement("COLLECTIONS")->LinkEndChild(xmlCol);
582 }
583 }
584 else // load file and read operator information
585 {
586 root = doc.FirstChildElement("NEKTAR");
587 xmlCol = root->FirstChildElement("COLLECTIONS");
589 xmlCol, m_timeLevel, false);
590 if (xmlCol)
591 {
592 // Read existing data
593 bool verbose = false;
594 ReadCollOps(xmlCol, global, verbose);
595 }
596 else
597 {
598 // Add timelevel tag
599 xmlCol = new TiXmlElement("TIMELEVEL");
600 xmlCol->SetAttribute("VALUE", m_timeLevel);
601 root->FirstChildElement("COLLECTIONS")->LinkEndChild(xmlCol);
602 }
603 }
604 }
605
606 // update global with m_opImpMap info
607 map<LibUtilities::ShapeType, int> ShapeMaxSize;
608 for (auto &opimp : m_opImpMap[m_timeLevel])
609 {
610 bool updateShape = true;
611 LibUtilities::ShapeType shape = opimp.first.GetShapeType();
612
613 // check to see if already added this shapes details but with
614 // a larger collection and if so do not update.
615 if (ShapeMaxSize.count(shape))
616 {
617 int ngeoms = opimp.first.GetNGeoms();
618 if (ngeoms > ShapeMaxSize[shape])
619 {
620 ShapeMaxSize[shape] = ngeoms;
621 }
622 else
623 {
624 updateShape = false;
625 }
626 }
627
628 if (updateShape)
629 {
630 for (auto &op : opimp.second)
631 {
632 global[op.first][ElmtOrder(shape, -1)] = op.second;
633 }
634 }
635 }
636
637 // loop over operators
638 for (auto &op : global)
639 {
640 // check to see which shapes are defined in this proc
641 Array<OneD, int> ElmtImp(LibUtilities::SIZE_ShapeType, -1);
642 Array<OneD, bool> ElmtDef(LibUtilities::SIZE_ShapeType, false);
643 for (auto &el : op.second)
644 {
645 ElmtImp[el.first.first] = el.second;
646 ElmtDef[el.first.first] = true;
647 }
648
649 comm->GetSpaceComm()->AllReduce(ElmtImp, LibUtilities::ReduceMax);
650
651 // loop over elements and update if not already defined
652 if (rank == 0)
653 {
654 for (int i = 1; i < LibUtilities::SIZE_ShapeType; ++i)
655 {
656 if ((ElmtImp[i] != -1) && (ElmtDef[i] == false))
657 {
658 global[op.first]
660 (ImplementationType)ElmtImp[i];
661 }
662 }
663 }
664 }
665
666 // Update Collection section with global data on root
667 if (rank == 0)
668 {
669 xmlCol->Clear();
670
671 map<LibUtilities::ShapeType, string> ShapeLetMap = {
679
680 for (auto &op : global)
681 {
682 TiXmlElement *ColOp = new TiXmlElement("OPERATOR");
683 xmlCol->LinkEndChild(ColOp);
684 ColOp->SetAttribute("TYPE", OperatorTypeMap[op.first]);
685
686 for (auto &el : op.second)
687 {
688 TiXmlElement *ElmtOp = new TiXmlElement("ELEMENT");
689 ColOp->LinkEndChild(ElmtOp);
690
691 ElmtOp->SetAttribute("TYPE", ShapeLetMap[el.first.first]);
692 ElmtOp->SetAttribute("ORDER", "*");
693 ElmtOp->SetAttribute("IMPTYPE",
694 ImplementationTypeMap[el.second]);
695 }
696 }
697
698 doc.SaveFile(outname);
699 }
700}
std::map< OperatorType, std::map< ElmtOrder, ImplementationType > > GlobalOpMap

References Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::ePyramid, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eSegment, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, Nektar::LibUtilities::SessionReader::GetXMLElementTimeLevel(), Nektar::Collections::ImplementationTypeMap, m_opImpMap, m_timeLevel, Nektar::Collections::OperatorTypeMap, ReadCollOps(), Nektar::LibUtilities::ReduceMax, and Nektar::LibUtilities::SIZE_ShapeType.

Referenced by Nektar::MultiRegions::ExpList::CreateCollections().

Member Data Documentation

◆ m_autotune

bool Nektar::Collections::CollectionOptimisation::m_autotune
private

Definition at line 168 of file CollectionOptimisation.h.

Referenced by CollectionOptimisation(), and IsUsingAutotuning().

◆ m_defaultType

ImplementationType Nektar::Collections::CollectionOptimisation::m_defaultType
private

◆ m_global

GlobalOpMap Nektar::Collections::CollectionOptimisation::m_global
private

Definition at line 167 of file CollectionOptimisation.h.

Referenced by CollectionOptimisation(), and GetOperatorImpMap().

◆ m_maxCollSize

size_t Nektar::Collections::CollectionOptimisation::m_maxCollSize
private

Definition at line 170 of file CollectionOptimisation.h.

Referenced by CollectionOptimisation(), and GetMaxCollectionSize().

◆ m_opImpMap

std::map< size_t, map< OpImpTimingKey, OperatorImpMap > > Nektar::Collections::CollectionOptimisation::m_opImpMap
staticprivate

Definition at line 166 of file CollectionOptimisation.h.

Referenced by SetWithTimings(), and UpdateOptFile().

◆ m_shapeDim

size_t Nektar::Collections::CollectionOptimisation::m_shapeDim
private

Definition at line 171 of file CollectionOptimisation.h.

Referenced by ReadCollOps().

◆ m_timeLevel

size_t Nektar::Collections::CollectionOptimisation::m_timeLevel
private

Definition at line 172 of file CollectionOptimisation.h.

Referenced by CollectionOptimisation(), SetWithTimings(), and UpdateOptFile().