Nektar++
Loading...
Searching...
No Matches
Movement/Movement.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: Movement.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: This file contains the base class for implementing
33// non-conformal geometry using the Movement object
34//
35////////////////////////////////////////////////////////////////////////////////
36
37#include <algorithm>
38#include <memory>
39#include <string>
40
45#include <tinyxml.h>
46
48{
49
50std::string static inline ReadTag(std::string &tagStr)
51{
52 std::string::size_type indxBeg = tagStr.find_first_of('[') + 1;
53 std::string::size_type indxEnd = tagStr.find_last_of(']') - 1;
54
56 indxBeg <= indxEnd,
57 (std::string("Error reading interface region definition:") + tagStr)
58 .c_str());
59
60 std::string indxStr = tagStr.substr(indxBeg, indxEnd - indxBeg + 1);
61
62 return indxStr;
63}
64
65std::string static inline StripParentheses(const std::string &str)
66{
67 auto length = str.length();
68 return str.substr(1, length - 2);
69}
70
72 MeshGraph *meshGraph)
73{
74 TiXmlNode *nektar = pSession->GetElement("NEKTAR");
75 if (nektar == nullptr)
76 {
77 return;
78 }
79
80 TiXmlNode *movement = nektar->FirstChild("MOVEMENT");
81 if (movement != nullptr)
82 {
83 bool zones = movement->FirstChild("ZONES") != nullptr;
84 if (zones)
85 {
86 TiXmlElement *zonesTag =
87 pSession->GetElement("NEKTAR/MOVEMENT/ZONES");
88 ReadZones(zonesTag, meshGraph, pSession);
89 }
90
91 bool interfaces = movement->FirstChild("INTERFACES") != nullptr;
92 if (interfaces)
93 {
94 TiXmlElement *interfacesTag =
95 pSession->GetElement("NEKTAR/MOVEMENT/INTERFACES");
96 ReadInterfaces(interfacesTag, meshGraph);
97 }
98
99 if (m_translate)
100 {
101 UpdateTransZoneBox(pSession);
102 }
103 }
104
105 // DEBUG COMMENTS
106 if (movement != nullptr && pSession->DefinesCmdLineArgument("verbose"))
107 {
108 auto comm = pSession->GetComm();
109 if (comm->TreatAsRankZero())
110 {
111 std::cout << "Movement Info:\n";
112 std::cout << "\tNum zones: " << m_zones.size() << "\n";
113 }
114
115 for (auto &zone : m_zones)
116 {
117 auto numEl = zone.second->GetElements().size();
118 comm->GetSpaceComm()->AllReduce(numEl, LibUtilities::ReduceSum);
119
120 // Find shape type if not on this proc
121 int shapeType =
122 !zone.second->GetElements().empty()
123 ? zone.second->GetElements().front()->GetShapeType()
124 : -1;
125 comm->GetSpaceComm()->AllReduce(shapeType, LibUtilities::ReduceMax);
126
127 if (comm->TreatAsRankZero())
128 {
129 std::cout << "\t- " << zone.first << " "
130 << MovementTypeStr[static_cast<int>(
131 zone.second->GetMovementType())]
132 << ": " << numEl << " "
133 << LibUtilities::ShapeTypeMap[shapeType] << "s\n";
134 }
135 }
136
137 if (comm->TreatAsRankZero())
138 {
139 std::cout << "\tNum interfaces: " << m_interfaces.size() << "\n";
140 }
141
142 for (auto &interface : m_interfaces)
143 {
144 auto numLeft =
145 interface.second->GetLeftInterface()->GetEdge().size();
146 auto numRight =
147 interface.second->GetRightInterface()->GetEdge().size();
148 comm->GetSpaceComm()->AllReduce(numLeft, LibUtilities::ReduceSum);
149 comm->GetSpaceComm()->AllReduce(numRight, LibUtilities::ReduceSum);
150
151 // Find shape type if not on this proc
152 int shapeTypeLeft =
153 !interface.second->GetLeftInterface()->GetEdge().empty()
154 ? interface.second->GetLeftInterface()
155 ->GetEdge()
156 .begin()
157 ->second->GetShapeType()
158 : -1;
159 comm->GetSpaceComm()->AllReduce(shapeTypeLeft,
161 int shapeTypeRight =
162 !interface.second->GetRightInterface()->GetEdge().empty()
163 ? interface.second->GetRightInterface()
164 ->GetEdge()
165 .begin()
166 ->second->GetShapeType()
167 : -1;
168 comm->GetSpaceComm()->AllReduce(shapeTypeRight,
170
171 if (comm->TreatAsRankZero())
172 {
173 std::cout << "\t- \"" << interface.first.second << "\": "
174 << interface.second->GetLeftInterface()->GetId()
175 << " (" << numLeft << " "
176 << LibUtilities::ShapeTypeMap[shapeTypeLeft]
177 << "s) <-> "
178 << interface.second->GetRightInterface()->GetId()
179 << " (" << numRight << " "
180 << LibUtilities::ShapeTypeMap[shapeTypeRight]
181 << "s)\n";
182 }
183 }
184
185 comm->GetSpaceComm()->Block();
186 if (comm->TreatAsRankZero())
187 {
188 std::cout << std::endl;
189 }
190 }
191}
192
193void Movement::ReadZones(TiXmlElement *zonesTag, MeshGraph *meshGraph,
195{
196 int coordDim = meshGraph->GetSpaceDimension();
197
198 ASSERTL0(zonesTag, "Unable to find ZONES tag in file.");
199 TiXmlElement *zonesElement = zonesTag->FirstChildElement();
200 while (zonesElement)
201 {
202 std::string zoneType = zonesElement->Value();
203
204 int err;
205 int indx;
206
207 err = zonesElement->QueryIntAttribute("ID", &indx);
208 ASSERTL0(err == TIXML_SUCCESS, "Unable to read zone ID.");
209
210 std::string interfaceDomainStr;
211 err = zonesElement->QueryStringAttribute("DOMAIN", &interfaceDomainStr);
212 ASSERTL0(err == TIXML_SUCCESS, "Unable to read zone domain.");
213
214 auto &domains = meshGraph->GetDomain();
215 auto domFind = std::stoi(ReadTag(interfaceDomainStr));
216 std::map<int, CompositeSharedPtr> domain;
217 if (domains.find(domFind) != domains.end())
218 {
219 domain = domains.at(domFind);
220 }
221
222 ZoneBaseShPtr zone;
223 if (zoneType == "F" || zoneType == "FIXED")
224 {
226 indx, domFind, domain, coordDim));
227 }
228 else if (zoneType == "R" || zoneType == "ROTATE" ||
229 zoneType == "ROTATING")
230 {
231 std::string originStr;
232 err = zonesElement->QueryStringAttribute("ORIGIN", &originStr);
233 ASSERTL0(err == TIXML_SUCCESS, "Unable to read origin.");
234 std::vector<NekDouble> originVec;
235 ParseUtils::GenerateVector(originStr, originVec);
236 NekPoint<NekDouble> origin =
237 NekPoint<NekDouble>(originVec[0], originVec[1], originVec[2]);
238
239 std::string axisStr;
240 err = zonesElement->QueryStringAttribute("AXIS", &axisStr);
241 ASSERTL0(err == TIXML_SUCCESS, "Unable to read axis.");
242 DNekVec axis(axisStr);
243 axis.Normalize();
244
245 std::string angularVelStr;
246 err = zonesElement->QueryStringAttribute("ANGVEL", &angularVelStr);
247 ASSERTL0(err == TIXML_SUCCESS, "Unable to read angular velocity.");
248
249 LibUtilities::EquationSharedPtr angularVelEqn =
251 pSession->GetInterpreter(), angularVelStr);
252
253 std::string rampTimeStr;
254 err = zonesElement->QueryStringAttribute("RAMPTIME", &rampTimeStr);
255 NekDouble rampTime = 0.0;
256 if (err == TIXML_SUCCESS)
257 {
258 LibUtilities::Interpreter expEvaluator;
259 int expr_id = expEvaluator.DefineFunction("", rampTimeStr);
260 rampTime = expEvaluator.Evaluate(expr_id);
261 }
262
263 std::string sectorStr;
264 err = zonesElement->QueryStringAttribute("SECTOR", &sectorStr);
265 NekDouble sector = 0.0;
267 if (err == TIXML_SUCCESS)
268 {
269
270 LibUtilities::Interpreter expEvaluator;
271 int expr_id = expEvaluator.DefineFunction("", sectorStr);
272 sector = expEvaluator.Evaluate(expr_id);
273
274 std::string baseStr;
275 err = zonesElement->QueryStringAttribute("BASE", &baseStr);
276 ASSERTL0(err == TIXML_SUCCESS, "Unable to read base.");
277 std::vector<NekDouble> baseVec;
278 ParseUtils::GenerateVector(baseStr, baseVec);
279 // Ensure the vector has exactly 3 elements before assignment
280 ASSERTL0(baseVec.size() == 3,
281 "BASE attribute must have exactly 3 elements.");
282 for (int i = 0; i < 3; ++i)
283 {
284 base[i] = baseVec[i];
285 }
286 m_sectorRotate = true;
287 }
288
290 indx, domFind, domain, coordDim, origin, axis, angularVelEqn,
291 rampTime, sector, base));
292
293 m_moveFlag = true;
294 m_rotate = true;
295 }
296 else if (zoneType == "T" || zoneType == "TRANSLATE" ||
297 zoneType == "TRANSLATING")
298 {
299 // Read velocity information
300 std::string velocityStr;
301 err = zonesElement->QueryStringAttribute("VELOCITY", &velocityStr);
302 ASSERTL0(err == TIXML_SUCCESS,
303 "Unable to read VELOCITY for translating zone " +
304 std::to_string(indx));
305
306 std::vector<std::string> velocityStrSplit;
307 ParseUtils::GenerateVector(velocityStr, velocityStrSplit);
308
309 ASSERTL0(velocityStrSplit.size() >= coordDim,
310 "VELOCITY dimensions must be greater than or equal to the "
311 "coordinate dimension for translating zone " +
312 std::to_string(indx));
313
315 for (int i = 0; i < coordDim; ++i)
316 {
317 velocityEqns[i] =
319 pSession->GetInterpreter(), velocityStrSplit[i]);
320 }
321
322 // Read displacement information
323 std::string displacementStr;
324 err = zonesElement->QueryStringAttribute("DISPLACEMENT",
325 &displacementStr);
326 ASSERTL0(err == TIXML_SUCCESS,
327 "Unable to read DISPLACEMENT for translating zone " +
328 std::to_string(indx));
329
330 std::vector<std::string> displacementStrSplit;
331 ParseUtils::GenerateVector(displacementStr, displacementStrSplit);
332
333 ASSERTL0(
334 displacementStrSplit.size() >= coordDim,
335 "DISPLACEMENT dimensions must be greater than or equal to the "
336 "coordinate dimension for translating zone " +
337 std::to_string(indx));
338
340 coordDim);
341 for (int i = 0; i < coordDim; ++i)
342 {
343 displacementEqns[i] =
345 pSession->GetInterpreter(), displacementStrSplit[i]);
346 }
347
348 zone = ZoneTranslateShPtr(
350 indx, domFind, domain, coordDim, velocityEqns,
351 displacementEqns));
352
353 m_moveFlag = true;
354 m_translate = true;
355 }
356 else
357 {
358 WARNINGL0(false, "Zone type '" + zoneType +
359 "' is unsupported. Valid types are: 'Fixed', "
360 "'Rotate', or 'Translate'.")
361 }
362
363 m_zones[indx] = zone;
364 zonesElement = zonesElement->NextSiblingElement();
365 }
366}
367
368void Movement::ReadInterfaces(TiXmlElement *interfacesTag, MeshGraph *meshGraph)
369{
370 ASSERTL0(interfacesTag, "Unable to find INTERFACES tag in file.");
371 TiXmlElement *interfaceElement = interfacesTag->FirstChildElement();
372
373 while (interfaceElement)
374 {
375 ASSERTL0(
376 "INTERFACE" == (std::string)interfaceElement->Value(),
377 "Only INTERFACE tags may be present inside the INTERFACES block.")
378
379 int err;
380
381 std::string name;
382 err = interfaceElement->QueryStringAttribute("NAME", &name);
383 ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface name.");
384 TiXmlElement *sideElement = interfaceElement->FirstChildElement();
385
386 Array<OneD, InterfaceShPtr> interfaces(2);
387 std::vector<int> cnt;
388 while (sideElement)
389 {
390 ASSERTL0(
391 cnt.size() < 2,
392 "In INTERFACE NAME " + name +
393 ", only two sides may be present in each INTERFACE block.")
394
395 // Read ID
396 int indx;
397 err = sideElement->QueryIntAttribute("ID", &indx);
398 ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface ID.");
399
400 // Read boundary
401 std::string boundaryStr;
402 err = sideElement->QueryStringAttribute("BOUNDARY", &boundaryStr);
403 ASSERTL0(err == TIXML_SUCCESS,
404 "Unable to read BOUNDARY in interface " +
405 std::to_string(indx));
406
407 CompositeMap boundaryEdge;
408 std::string indxStr = ReadTag(boundaryStr);
409 meshGraph->GetCompositeList(indxStr, boundaryEdge);
410
411 // Read SKIPCHECK flag
412 bool skipCheck = false;
413 err = sideElement->QueryBoolAttribute("SKIPCHECK", &skipCheck);
414
415 // Sets location in interface pair to 0 for left, and 1 for right
416 auto sideElVal = sideElement->ValueStr();
417 if (sideElVal == "LEFT" || sideElVal == "L")
418 {
419 cnt.emplace_back(0);
420 }
421 else if (sideElVal == "RIGHT" || sideElVal == "R")
422 {
423 cnt.emplace_back(1);
424 }
425 else
426 {
428 sideElement->ValueStr() +
429 " is not a valid interface side for interface "
430 "NAME " +
431 name + ". Please only use LEFT or RIGHT.")
432 }
433
434 interfaces[cnt[cnt.size() - 1]] =
436 indx, boundaryEdge, skipCheck));
437
438 sideElement = sideElement->NextSiblingElement();
439 }
440
441 ASSERTL0(std::accumulate(cnt.begin(), cnt.end(), 0) == 1,
442 "You must have only one LEFT and one RIGHT side"
443 " present in interface NAME " +
444 name)
445
446 m_interfaces[std::make_pair(m_interfaces.size(), name)] =
448 interfaces[0], interfaces[1]));
449 interfaceElement = interfaceElement->NextSiblingElement();
450 }
451}
452
453/// Export this Movement information to a Nektar++ XML file.
454void Movement::WriteMovement(TiXmlElement *root)
455{
456 if (m_zones.size() == 0 && m_interfaces.size() == 0)
457 {
458 return;
459 }
460 TiXmlElement *movement = new TiXmlElement("MOVEMENT");
461 root->LinkEndChild(movement);
462
463 TiXmlElement *zones = new TiXmlElement("ZONES");
464 for (auto &i : m_zones)
465 {
466 const ZoneBaseShPtr z = i.second;
467 const MovementType mtype = z->GetMovementType();
468 std::string label = MovementTypeStr[static_cast<int>(mtype)];
469 TiXmlElement *e = new TiXmlElement(label.substr(0, 1));
470 e->SetAttribute("ID", i.first);
471 std::stringstream s;
472 s << "D[" << z->GetDomainID() << "]";
473 e->SetAttribute("DOMAIN", s.str());
474
475 switch (mtype)
476 {
478 {
479 auto rotate = std::static_pointer_cast<ZoneRotate>(z);
480 e->SetAttribute(
481 "ORIGIN", StripParentheses(rotate->GetOrigin().AsString()));
482 e->SetAttribute("AXIS",
483 StripParentheses(rotate->GetAxis().AsString()));
484 e->SetAttribute("ANGVEL",
485 rotate->GetAngularVelEqn()->GetExpression());
486 }
487 break;
489 {
490 auto translate = std::static_pointer_cast<ZoneTranslate>(z);
491 e->SetAttribute(
492 "XVELOCITY",
493 translate->GetVelocityEquation()[0]->GetExpression());
494 e->SetAttribute(
495 "XDISPLACEMENT",
496 translate->GetDisplacementEquation()[0]->GetExpression());
497 e->SetAttribute(
498 "YVELOCITY",
499 translate->GetVelocityEquation()[1]->GetExpression());
500 e->SetAttribute(
501 "YDISPLACEMENT",
502 translate->GetDisplacementEquation()[1]->GetExpression());
503 e->SetAttribute(
504 "ZVELOCITY",
505 translate->GetVelocityEquation()[2]->GetExpression());
506 e->SetAttribute(
507 "ZDISPLACEMENT",
508 translate->GetDisplacementEquation()[2]->GetExpression());
509 }
510 break;
511 default:
512 break;
513 }
514 zones->LinkEndChild(e);
515 }
516 movement->LinkEndChild(zones);
517
518 TiXmlElement *interfaces = new TiXmlElement("INTERFACES");
519 for (auto &i : m_interfaces)
520 {
521 const std::string interfaceName = i.first.second;
522 TiXmlElement *e = new TiXmlElement("INTERFACE");
523 e->SetAttribute("NAME", interfaceName);
524 const InterfaceShPtr left = i.second->GetLeftInterface();
525 const InterfaceShPtr right = i.second->GetRightInterface();
526 if (left)
527 {
528 TiXmlElement *left_e = new TiXmlElement("L");
529 left_e->SetAttribute("ID", left->GetId());
530 left_e->SetAttribute(
531 "BOUNDARY",
532 "C[" + ParseUtils::GenerateSeqString(left->GetCompositeIDs()) +
533 "]");
534 e->LinkEndChild(left_e);
535 }
536 if (right)
537 {
538 TiXmlElement *right_e = new TiXmlElement("R");
539 right_e->SetAttribute("ID", right->GetId());
540 right_e->SetAttribute(
541 "BOUNDARY",
542 "C[" + ParseUtils::GenerateSeqString(right->GetCompositeIDs()) +
543 "]");
544 e->LinkEndChild(right_e);
545 }
546 interfaces->LinkEndChild(e);
547 }
548 movement->LinkEndChild(interfaces);
549}
550
551// Acts as a placeholder for when ALE function and moving geometry capability
552// is added. Currently unused.
554{
555 std::set<int> movedZoneIds;
556 for (auto &zone : m_zones)
557 {
558 if (zone.second->Move(timeStep))
559 {
560 movedZoneIds.insert(zone.first);
561 }
562 }
563
564 // If zone has moved, set all interfaces on that zone to moved.
565 for (auto &interPair : m_interfaces)
566 {
567 int leftId = interPair.second->GetLeftInterface()->GetId();
568 int rightId = interPair.second->GetRightInterface()->GetId();
569
570 if (movedZoneIds.find(leftId) != movedZoneIds.end() ||
571 movedZoneIds.find(rightId) != movedZoneIds.end())
572 {
573 m_zones[leftId]->GetMoved() = true;
574 m_zones[rightId]->GetMoved() = true;
575 }
576 }
577 m_moved = true;
578}
579
580/// Store a zone object with this Movement data
582{
583 m_zones[zone->GetId()] = zone;
584 MovementType mtype = zone->GetMovementType();
585 if (mtype != MovementType::eFixed && mtype != MovementType::eNone)
586 {
587 m_moveFlag = true;
588 }
589}
590
591/// Store an interface pair with this Movement data
592void Movement::AddInterface(std::string name, InterfaceShPtr left,
593 InterfaceShPtr right)
594{
595 m_interfaces[std::make_pair(m_interfaces.size(), name)] =
598}
599
600/// Generate domain box for translation mesh
603{
604 // Loop over all translate zones to get the box of the whole zone
605 for (auto &zones : m_zones)
606 {
607 if (zones.second->GetMovementType() == eTranslate)
608 {
609 auto zone = std::static_pointer_cast<ZoneTranslate>(zones.second);
610 Array<OneD, NekDouble> ZoneBox(6);
611 Array<OneD, NekDouble> ZoneLength(3);
612 ZoneBox = zone->GetZoneBox();
613 ZoneLength = zone->GetZoneLength();
614 auto comm = pSession->GetComm();
615 for (int i = 0; i < 3; ++i)
616 {
617 auto &min = ZoneBox[i];
618 auto &max = ZoneBox[i + 3];
619 comm->GetSpaceComm()->AllReduce(min, LibUtilities::ReduceMin);
620 comm->GetSpaceComm()->AllReduce(max, LibUtilities::ReduceMax);
621 ZoneLength[i] = max - min;
622 }
623 // Update the zone box for all translate zones
624 zone->UpdateZoneBox(ZoneBox, ZoneLength);
625 }
626 }
627}
628
629} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define WARNINGL0(condition, msg)
Interpreter class for the evaluation of mathematical expressions.
Definition Interpreter.h:76
int DefineFunction(const std::string &vlist, const std::string &expr)
Defines a function for the purposes of evaluation.
NekDouble Evaluate(const int id)
Evaluate a function which depends only on constants and/or parameters.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static std::string GenerateSeqString(const std::vector< T > &v)
Generate a compressed comma-separated string representation of a vector of unsigned integers.
Definition ParseUtils.h:72
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Base class for a spectral/hp element mesh.
Definition MeshGraph.h:290
std::map< int, std::map< int, CompositeSharedPtr > > & GetDomain()
Definition MeshGraph.h:383
void GetCompositeList(const std::string &compositeStr, CompositeMap &compositeVector) const
int GetSpaceDimension()
Dimension of the space (can be a 1D curve in 3D space).
Definition MeshGraph.h:323
void UpdateTransZoneBox(const LibUtilities::SessionReaderSharedPtr &pSession)
Update the box of translate zones.
void PerformMovement(NekDouble timeStep)
void ReadZones(TiXmlElement *zonesTag, MeshGraph *meshGraph, const LibUtilities::SessionReaderSharedPtr &pSession)
Read zones given TiXmlDocument.
void AddInterface(std::string name, InterfaceShPtr left, InterfaceShPtr right)
Add pair of interfaces to this data.
void ReadInterfaces(TiXmlElement *interfacesTag, MeshGraph *meshGraph)
Read interfaces given TiXmlDocument.
void AddZone(ZoneBaseShPtr zone)
Add a zone object to this Movement data.
void WriteMovement(TiXmlElement *root)
Write the MOVEMENT section of the XML file.
InterfaceCollection m_interfaces
Definition Movement.h:141
std::map< int, ZoneBaseShPtr > m_zones
Definition Movement.h:142
Movement()=default
Default constructor.
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition ShapeType.hpp:75
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition Equation.h:131
static std::string StripParentheses(const std::string &str)
std::shared_ptr< ZoneBase > ZoneBaseShPtr
Definition Zones.h:177
std::shared_ptr< ZoneTranslate > ZoneTranslateShPtr
Definition Zones.h:376
std::shared_ptr< InterfacePair > InterfacePairShPtr
std::shared_ptr< ZoneRotate > ZoneRotateShPtr
Definition Zones.h:375
std::shared_ptr< Interface > InterfaceShPtr
MovementType
Enum of zone movement type.
Definition Zones.h:48
std::shared_ptr< ZoneFixed > ZoneFixedShPtr
Definition Zones.h:377
static std::string ReadTag(std::string &tagStr)
const std::string MovementTypeStr[]
Map of zone movement type to movement type string.
Definition Zones.h:57
std::map< int, CompositeSharedPtr > CompositeMap
Definition MeshGraph.h:179
scalarT< T > max(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:305
scalarT< T > min(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:300