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 ()=default
 
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 119 of file CollectionOptimisation.h.

Member Typedef Documentation

◆ ElmtOrder

Definition at line 158 of file CollectionOptimisation.h.

◆ GlobalOpMap

Definition at line 160 of file CollectionOptimisation.h.

Constructor & Destructor Documentation

◆ CollectionOptimisation()

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

Definition at line 55 of file CollectionOptimisation.cpp.

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

References ASSERTL0, Nektar::Collections::eHelmholtz, Nektar::LibUtilities::eHexahedron, Nektar::Collections::eIterPerExp, Nektar::Collections::eMatrixFree, Nektar::Collections::eNoCollection, Nektar::Collections::eNoImpType, Nektar::Collections::ePhysDeriv, Nektar::Collections::ePhysInterp1DScaled, 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 ( )
default

Member Function Documentation

◆ GetDefaultImplementationType()

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

Definition at line 129 of file CollectionOptimisation.h.

130 {
131 return m_defaultType;
132 }

References m_defaultType.

◆ GetMaxCollectionSize()

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

Definition at line 134 of file CollectionOptimisation.h.

135 {
136 return m_maxCollSize;
137 }

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 356 of file CollectionOptimisation.cpp.

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

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 139 of file CollectionOptimisation.h.

140 {
141 return m_autotune;
142 }

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 227 of file CollectionOptimisation.cpp.

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

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 393 of file CollectionOptimisation.cpp.

396{
397 OperatorImpMap ret;
398
399 StdRegions::StdExpansionSharedPtr pExp = pCollExp[0];
400
401 // check to see if already defined for this expansion
402 OpImpTimingKey OpKey(pExp, pCollExp.size(), pExp->GetNumBases());
403 if (m_opImpMap.count(m_timeLevel) != 0 &&
404 m_opImpMap[m_timeLevel].count(OpKey) != 0)
405 {
406 ret = m_opImpMap[m_timeLevel][OpKey];
407 return ret;
408 }
409
410 LibUtilities::Timer t;
411
412 if (verbose)
413 {
414 cout << "Collection Implementation for "
415 << LibUtilities::ShapeTypeMap[pExp->DetShapeType()] << " ( ";
416 for (int i = 0; i < pExp->GetNumBases(); ++i)
417 {
418 cout << pExp->GetBasis(i)->GetNumModes() << " ";
419 }
420 cout << ")"
421 << " for ngeoms = " << pCollExp.size() << endl;
422 }
423 // set up an array of collections
424 CollectionVector coll;
425
426 StdRegions::ConstFactorMap factors; // required for helmholtz operator
428
429 // do we need to dynamically update the value of this eFactorconst based on
430 // the solver that is called?
432 int qpInsideElmt_idir; // declaration for the total number of points towards
433 // the i-th direction
434 int newQpInsideElmt_idir; // declaration for the scaled total number of
435 // points towards the i-th direction
436 int newTotQpInsideElmt{1}; // initialization for the total number of scaled
437 // quadrature points inside an element
438 for (int i = 0; i < pExp->GetShapeDimension(); ++i)
439 {
440 qpInsideElmt_idir = pExp->GetNumPoints(i);
441 newQpInsideElmt_idir =
442 (int)qpInsideElmt_idir * factors[StdRegions::eFactorConst];
443 newTotQpInsideElmt *= newQpInsideElmt_idir;
444 }
445 int maxsize = pCollExp.size() * max(pExp->GetNcoeffs(), newTotQpInsideElmt);
446 Array<OneD, NekDouble> inarray(maxsize, 1.0);
447 Array<OneD, NekDouble> outarray1(maxsize);
448 Array<OneD, NekDouble> outarray2(maxsize);
449 Array<OneD, NekDouble> outarray3(maxsize);
450
451 // Advection velocities required for optimisation of linearADR operator
452 StdRegions::VarCoeffMap varcoeffs;
456 for (int i = 0; i < pExp->GetShapeDimension(); i++)
457 {
458 varcoeffs[varcoefftypes[i]] = Array<OneD, NekDouble>(maxsize, 1.0);
459 }
460
461 for (int imp = 1; imp < SIZE_ImplementationType; ++imp)
462 {
464 OperatorImpMap impTypes;
465 std::vector<OperatorType> opTypes;
466 for (int i = 0; i < SIZE_OperatorType; ++i)
467 {
468 OperatorType opType = (OperatorType)i;
469 OperatorKey opKey(pCollExp[0]->DetShapeType(), opType, impType,
470 pCollExp[0]->IsNodalNonTensorialExp());
471
472 if (GetOperatorFactory().ModuleExists(opKey))
473 {
474 impTypes[opType] = impType;
475 opTypes.push_back(opType);
476 }
477 }
478
479 Collection collLoc(pCollExp, impTypes);
480 for (auto opType : opTypes)
481 {
482 collLoc.Initialise(opType, factors);
483 // Add varcoeffs only for ADR operator
485 {
486 collLoc.UpdateVarcoeffs(opType, varcoeffs);
487 }
488 }
489 coll.push_back(collLoc);
490 }
491
492 // Determine the number of tests to do in one second
493 Array<OneD, int> Ntest(SIZE_OperatorType);
494 for (int i = 0; i < SIZE_OperatorType; ++i)
495 {
496 OperatorType OpType = (OperatorType)i;
497
498 t.Start();
499
500 coll[0].ApplyOperator(OpType, inarray, outarray1, outarray2, outarray3);
501 t.Stop();
502
503 NekDouble oneTest = t.TimePerTest(1);
504
505 Ntest[i] = max((int)(0.25 / oneTest), 1);
506 }
507
508 Array<OneD, NekDouble> timing(SIZE_ImplementationType);
509
510 if (verbose)
511 {
512 cout << "\t "
513 << " Op. "
514 << ":\t"
515 << "opt. Impl."
516 << "\t (IterLocExp, IterStdExp, "
517 "StdMat, SumFac, MatrixFree)"
518 << endl;
519 }
520
521 // loop over all operators and determine fastest implementation
522 for (int i = 0; i < SIZE_OperatorType; ++i)
523 {
524 OperatorType OpType = (OperatorType)i;
525
526 // call collection implementation in thorugh ExpList.
527 for (int imp = 0; imp < coll.size(); ++imp)
528 {
529 if (coll[imp].HasOperator(OpType))
530 {
531 t.Start();
532 for (int n = 0; n < Ntest[i]; ++n)
533 {
534 coll[imp].ApplyOperator(OpType, inarray, outarray1,
535 outarray2, outarray3);
536 }
537 t.Stop();
538 timing[imp] = t.TimePerTest(Ntest[i]);
539 }
540 else
541 {
542 timing[imp] = 1000.0;
543 }
544 }
545 // determine optimal implementation. Note +1 to
546 // remove NoImplementationType flag
547 int minImp = Vmath::Imin(coll.size(), timing, 1) + 1;
548
549 if (verbose)
550 {
551 cout << "\t " << OperatorTypeMap1[i] << ": \t"
552 << ImplementationTypeMap1[minImp] << "\t (";
553 for (int j = 0; j < coll.size(); ++j)
554 {
555 if (timing[j] > 999.0)
556 {
557 cout << " -- ";
558 }
559 else
560 {
561 cout << timing[j];
562 }
563 if (j != coll.size() - 1)
564 {
565 cout << ", ";
566 }
567 }
568 cout << ")" << endl;
569 }
570
571 // set up new map
572 ret[OpType] = (ImplementationType)minImp;
573 }
574
575 // store map for use by another expansion.
576 if (m_opImpMap.count(m_timeLevel) == 0)
577 {
578 m_opImpMap[m_timeLevel] = map<OpImpTimingKey, OperatorImpMap>();
579 }
580 m_opImpMap[m_timeLevel][OpKey] = ret;
581 return ret;
582}
static std::map< size_t, std::map< OpImpTimingKey, OperatorImpMap > > m_opImpMap
const char *const ImplementationTypeMap1[]
Definition: Operator.h:104
@ eLinearAdvectionDiffusionReaction
Definition: Operator.h:66
std::tuple< LibUtilities::ShapeType, OperatorType, ImplementationType, ExpansionIsNodal > OperatorKey
Key for describing an Operator.
Definition: Operator.h:120
OperatorFactory & GetOperatorFactory()
Returns the singleton Operator factory object.
Definition: Operator.cpp:44
const char *const OperatorTypeMap1[]
Definition: Operator.h:82
std::vector< Collection > CollectionVector
Definition: Collection.h:128
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375
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.hpp:704

References Nektar::StdRegions::eFactorConst, Nektar::StdRegions::eFactorLambda, Nektar::Collections::eLinearAdvectionDiffusionReaction, Nektar::StdRegions::eVarCoeffVelX, Nektar::StdRegions::eVarCoeffVelY, Nektar::StdRegions::eVarCoeffVelZ, 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(), Nektar::LibUtilities::Timer::TimePerTest(), and Nektar::Collections::Collection::UpdateVarcoeffs().

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

◆ UpdateOptFile()

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

Definition at line 584 of file CollectionOptimisation.cpp.

586{
587 if (comm->IsParallelInTime() && comm->GetTimeComm()->GetRank() > 0)
588 {
589 // No need to repeatly update the optfile for each time chunk.
590 return;
591 }
592
593 std::string outname = sessName.substr(0, sessName.find("_xml/")) + ".opt";
594
595 TiXmlDocument doc;
596 TiXmlElement *root;
597 TiXmlElement *xmlCol = new TiXmlElement("COLLECTIONS");
598 GlobalOpMap global;
599 int rank = comm->GetSpaceComm()->GetRank();
600 if (rank == 0)
601 {
602 if (!doc.LoadFile(outname)) // set up new file
603 {
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())
610 {
611 // Add timelevel tag
612 xmlCol = new TiXmlElement("TIMELEVEL");
613 xmlCol->SetAttribute("VALUE", 0);
614 root->FirstChildElement("COLLECTIONS")->LinkEndChild(xmlCol);
615 }
616 }
617 else // load file and read operator information
618 {
619 root = doc.FirstChildElement("NEKTAR");
620 xmlCol = root->FirstChildElement("COLLECTIONS");
622 xmlCol, m_timeLevel, false);
623 if (xmlCol)
624 {
625 // Read existing data
626 bool verbose = false;
627 ReadCollOps(xmlCol, global, verbose);
628 }
629 else
630 {
631 // Add timelevel tag
632 xmlCol = new TiXmlElement("TIMELEVEL");
633 xmlCol->SetAttribute("VALUE", m_timeLevel);
634 root->FirstChildElement("COLLECTIONS")->LinkEndChild(xmlCol);
635 }
636 }
637 }
638
639 // update global with m_opImpMap info
640 map<LibUtilities::ShapeType, int> ShapeMaxSize;
641 for (auto &opimp : m_opImpMap[m_timeLevel])
642 {
643 bool updateShape = true;
644 LibUtilities::ShapeType shape = opimp.first.GetShapeType();
645
646 // check to see if already added this shapes details but with
647 // a larger collection and if so do not update.
648 if (ShapeMaxSize.count(shape))
649 {
650 int ngeoms = opimp.first.GetNGeoms();
651 if (ngeoms > ShapeMaxSize[shape])
652 {
653 ShapeMaxSize[shape] = ngeoms;
654 }
655 else
656 {
657 updateShape = false;
658 }
659 }
660
661 if (updateShape)
662 {
663 for (auto &op : opimp.second)
664 {
665 global[op.first][ElmtOrder(shape, -1)] = op.second;
666 }
667 }
668 }
669
670 // loop over operators
671 for (auto &op : global)
672 {
673 // check to see which shapes are defined in this proc
674 Array<OneD, int> ElmtImp(LibUtilities::SIZE_ShapeType, -1);
675 Array<OneD, bool> ElmtDef(LibUtilities::SIZE_ShapeType, false);
676 for (auto &el : op.second)
677 {
678 ElmtImp[el.first.first] = el.second;
679 ElmtDef[el.first.first] = true;
680 }
681
682 comm->GetSpaceComm()->AllReduce(ElmtImp, LibUtilities::ReduceMax);
683
684 // loop over elements and update if not already defined
685 if (rank == 0)
686 {
687 for (int i = 1; i < LibUtilities::SIZE_ShapeType; ++i)
688 {
689 if ((ElmtImp[i] != -1) && (ElmtDef[i] == false))
690 {
691 global[op.first]
693 (ImplementationType)ElmtImp[i];
694 }
695 }
696 }
697 }
698
699 // Update Collection section with global data on root
700 if (rank == 0)
701 {
702 xmlCol->Clear();
703
704 map<LibUtilities::ShapeType, string> ShapeLetMap = {
712
713 for (auto &op : global)
714 {
715 TiXmlElement *ColOp = new TiXmlElement("OPERATOR");
716 xmlCol->LinkEndChild(ColOp);
717 ColOp->SetAttribute("TYPE", OperatorTypeMap[op.first]);
718
719 for (auto &el : op.second)
720 {
721 TiXmlElement *ElmtOp = new TiXmlElement("ELEMENT");
722 ColOp->LinkEndChild(ElmtOp);
723
724 ElmtOp->SetAttribute("TYPE", ShapeLetMap[el.first.first]);
725 ElmtOp->SetAttribute("ORDER", "*");
726 ElmtOp->SetAttribute("IMPTYPE",
727 ImplementationTypeMap[el.second]);
728 }
729 }
730
731 doc.SaveFile(outname);
732 }
733}
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 164 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 163 of file CollectionOptimisation.h.

Referenced by CollectionOptimisation(), and GetOperatorImpMap().

◆ m_maxCollSize

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

Definition at line 166 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 162 of file CollectionOptimisation.h.

Referenced by SetWithTimings(), and UpdateOptFile().

◆ m_shapeDim

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

Definition at line 167 of file CollectionOptimisation.h.

Referenced by ReadCollOps().

◆ m_timeLevel

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

Definition at line 168 of file CollectionOptimisation.h.

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