50std::string
static inline ReadTag(std::string &tagStr)
52 std::string::size_type indxBeg = tagStr.find_first_of(
'[') + 1;
53 std::string::size_type indxEnd = tagStr.find_last_of(
']') - 1;
57 (std::string(
"Error reading interface region definition:") + tagStr)
60 std::string indxStr = tagStr.substr(indxBeg, indxEnd - indxBeg + 1);
67 auto length = str.length();
68 return str.substr(1, length - 2);
74 TiXmlNode *
nektar = pSession->GetElement(
"NEKTAR");
80 TiXmlNode *movement =
nektar->FirstChild(
"MOVEMENT");
81 if (movement !=
nullptr)
83 bool zones = movement->FirstChild(
"ZONES") !=
nullptr;
86 TiXmlElement *zonesTag =
87 pSession->GetElement(
"NEKTAR/MOVEMENT/ZONES");
91 bool interfaces = movement->FirstChild(
"INTERFACES") !=
nullptr;
94 TiXmlElement *interfacesTag =
95 pSession->GetElement(
"NEKTAR/MOVEMENT/INTERFACES");
106 if (movement !=
nullptr && pSession->DefinesCmdLineArgument(
"verbose"))
108 auto comm = pSession->GetComm();
109 if (comm->TreatAsRankZero())
111 std::cout <<
"Movement Info:\n";
112 std::cout <<
"\tNum zones: " <<
m_zones.size() <<
"\n";
117 auto numEl = zone.second->GetElements().size();
122 !zone.second->GetElements().empty()
123 ? zone.second->GetElements().front()->GetShapeType()
127 if (comm->TreatAsRankZero())
129 std::cout <<
"\t- " << zone.first <<
" "
131 zone.second->GetMovementType())]
132 <<
": " << numEl <<
" "
137 if (comm->TreatAsRankZero())
139 std::cout <<
"\tNum interfaces: " <<
m_interfaces.size() <<
"\n";
145 interface.second->GetLeftInterface()->GetEdge().size();
147 interface.second->GetRightInterface()->GetEdge().size();
153 !interface.second->GetLeftInterface()->GetEdge().empty()
154 ? interface.second->GetLeftInterface()
157 ->second->GetShapeType()
159 comm->GetSpaceComm()->AllReduce(shapeTypeLeft,
162 !interface.second->GetRightInterface()->GetEdge().empty()
163 ? interface.second->GetRightInterface()
166 ->second->GetShapeType()
168 comm->GetSpaceComm()->AllReduce(shapeTypeRight,
171 if (comm->TreatAsRankZero())
173 std::cout <<
"\t- \"" << interface.first.second <<
"\": "
174 << interface.second->GetLeftInterface()->GetId()
175 <<
" (" << numLeft <<
" "
178 << interface.second->GetRightInterface()->GetId()
179 <<
" (" << numRight <<
" "
185 comm->GetSpaceComm()->Block();
186 if (comm->TreatAsRankZero())
188 std::cout << std::endl;
198 ASSERTL0(zonesTag,
"Unable to find ZONES tag in file.");
199 TiXmlElement *zonesElement = zonesTag->FirstChildElement();
202 std::string zoneType = zonesElement->Value();
207 err = zonesElement->QueryIntAttribute(
"ID", &indx);
208 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read zone ID.");
210 std::string interfaceDomainStr;
211 err = zonesElement->QueryStringAttribute(
"DOMAIN", &interfaceDomainStr);
212 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read zone domain.");
215 auto domFind = std::stoi(
ReadTag(interfaceDomainStr));
216 std::map<int, CompositeSharedPtr> domain;
217 if (domains.find(domFind) != domains.end())
219 domain = domains.at(domFind);
223 if (zoneType ==
"F" || zoneType ==
"FIXED")
226 indx, domFind, domain, coordDim));
228 else if (zoneType ==
"R" || zoneType ==
"ROTATE" ||
229 zoneType ==
"ROTATING")
231 std::string originStr;
232 err = zonesElement->QueryStringAttribute(
"ORIGIN", &originStr);
233 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read origin.");
234 std::vector<NekDouble> originVec;
240 err = zonesElement->QueryStringAttribute(
"AXIS", &axisStr);
241 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read axis.");
245 std::string angularVelStr;
246 err = zonesElement->QueryStringAttribute(
"ANGVEL", &angularVelStr);
247 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read angular velocity.");
251 pSession->GetInterpreter(), angularVelStr);
253 std::string rampTimeStr;
254 err = zonesElement->QueryStringAttribute(
"RAMPTIME", &rampTimeStr);
256 if (err == TIXML_SUCCESS)
260 rampTime = expEvaluator.
Evaluate(expr_id);
263 std::string sectorStr;
264 err = zonesElement->QueryStringAttribute(
"SECTOR", §orStr);
267 if (err == TIXML_SUCCESS)
272 sector = expEvaluator.
Evaluate(expr_id);
275 err = zonesElement->QueryStringAttribute(
"BASE", &baseStr);
276 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read base.");
277 std::vector<NekDouble> baseVec;
281 "BASE attribute must have exactly 3 elements.");
282 for (
int i = 0; i < 3; ++i)
284 base[i] = baseVec[i];
290 indx, domFind, domain, coordDim, origin, axis, angularVelEqn,
291 rampTime, sector, base));
296 else if (zoneType ==
"T" || zoneType ==
"TRANSLATE" ||
297 zoneType ==
"TRANSLATING")
300 std::string velocityStr;
301 err = zonesElement->QueryStringAttribute(
"VELOCITY", &velocityStr);
303 "Unable to read VELOCITY for translating zone " +
304 std::to_string(indx));
306 std::vector<std::string> velocityStrSplit;
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));
315 for (
int i = 0; i < coordDim; ++i)
319 pSession->GetInterpreter(), velocityStrSplit[i]);
323 std::string displacementStr;
324 err = zonesElement->QueryStringAttribute(
"DISPLACEMENT",
327 "Unable to read DISPLACEMENT for translating zone " +
328 std::to_string(indx));
330 std::vector<std::string> displacementStrSplit;
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));
341 for (
int i = 0; i < coordDim; ++i)
343 displacementEqns[i] =
345 pSession->GetInterpreter(), displacementStrSplit[i]);
350 indx, domFind, domain, coordDim, velocityEqns,
358 WARNINGL0(
false,
"Zone type '" + zoneType +
359 "' is unsupported. Valid types are: 'Fixed', "
360 "'Rotate', or 'Translate'.")
364 zonesElement = zonesElement->NextSiblingElement();
370 ASSERTL0(interfacesTag,
"Unable to find INTERFACES tag in file.");
371 TiXmlElement *interfaceElement = interfacesTag->FirstChildElement();
373 while (interfaceElement)
376 "INTERFACE" == (std::string)interfaceElement->Value(),
377 "Only INTERFACE tags may be present inside the INTERFACES block.")
382 err = interfaceElement->QueryStringAttribute(
"NAME", &name);
383 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read interface name.");
384 TiXmlElement *sideElement = interfaceElement->FirstChildElement();
387 std::vector<int> cnt;
392 "In INTERFACE NAME " + name +
393 ", only two sides may be present in each INTERFACE block.")
397 err = sideElement->QueryIntAttribute(
"ID", &indx);
398 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read interface ID.");
401 std::string boundaryStr;
402 err = sideElement->QueryStringAttribute(
"BOUNDARY", &boundaryStr);
404 "Unable to read BOUNDARY in interface " +
405 std::to_string(indx));
408 std::string indxStr =
ReadTag(boundaryStr);
412 bool skipCheck =
false;
413 err = sideElement->QueryBoolAttribute(
"SKIPCHECK", &skipCheck);
416 auto sideElVal = sideElement->ValueStr();
417 if (sideElVal ==
"LEFT" || sideElVal ==
"L")
421 else if (sideElVal ==
"RIGHT" || sideElVal ==
"R")
428 sideElement->ValueStr() +
429 " is not a valid interface side for interface "
431 name +
". Please only use LEFT or RIGHT.")
434 interfaces[cnt[cnt.size() - 1]] =
436 indx, boundaryEdge, skipCheck));
438 sideElement = sideElement->NextSiblingElement();
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 " +
448 interfaces[0], interfaces[1]));
449 interfaceElement = interfaceElement->NextSiblingElement();
460 TiXmlElement *movement =
new TiXmlElement(
"MOVEMENT");
461 root->LinkEndChild(movement);
463 TiXmlElement *zones =
new TiXmlElement(
"ZONES");
469 TiXmlElement *e =
new TiXmlElement(label.substr(0, 1));
470 e->SetAttribute(
"ID", i.first);
472 s <<
"D[" << z->GetDomainID() <<
"]";
473 e->SetAttribute(
"DOMAIN", s.str());
479 auto rotate = std::static_pointer_cast<ZoneRotate>(z);
482 e->SetAttribute(
"AXIS",
484 e->SetAttribute(
"ANGVEL",
485 rotate->GetAngularVelEqn()->GetExpression());
490 auto translate = std::static_pointer_cast<ZoneTranslate>(z);
493 translate->GetVelocityEquation()[0]->GetExpression());
496 translate->GetDisplacementEquation()[0]->GetExpression());
499 translate->GetVelocityEquation()[1]->GetExpression());
502 translate->GetDisplacementEquation()[1]->GetExpression());
505 translate->GetVelocityEquation()[2]->GetExpression());
508 translate->GetDisplacementEquation()[2]->GetExpression());
514 zones->LinkEndChild(e);
516 movement->LinkEndChild(zones);
518 TiXmlElement *interfaces =
new TiXmlElement(
"INTERFACES");
521 const std::string interfaceName = i.first.second;
522 TiXmlElement *e =
new TiXmlElement(
"INTERFACE");
523 e->SetAttribute(
"NAME", interfaceName);
528 TiXmlElement *left_e =
new TiXmlElement(
"L");
529 left_e->SetAttribute(
"ID", left->GetId());
530 left_e->SetAttribute(
534 e->LinkEndChild(left_e);
538 TiXmlElement *right_e =
new TiXmlElement(
"R");
539 right_e->SetAttribute(
"ID", right->GetId());
540 right_e->SetAttribute(
544 e->LinkEndChild(right_e);
546 interfaces->LinkEndChild(e);
548 movement->LinkEndChild(interfaces);
555 std::set<int> movedZoneIds;
558 if (zone.second->Move(timeStep))
560 movedZoneIds.insert(zone.first);
567 int leftId = interPair.second->GetLeftInterface()->GetId();
568 int rightId = interPair.second->GetRightInterface()->GetId();
570 if (movedZoneIds.find(leftId) != movedZoneIds.end() ||
571 movedZoneIds.find(rightId) != movedZoneIds.end())
573 m_zones[leftId]->GetMoved() =
true;
574 m_zones[rightId]->GetMoved() =
true;
607 if (zones.second->GetMovementType() ==
eTranslate)
609 auto zone = std::static_pointer_cast<ZoneTranslate>(zones.second);
612 ZoneBox = zone->GetZoneBox();
613 ZoneLength = zone->GetZoneLength();
614 auto comm = pSession->GetComm();
615 for (
int i = 0; i < 3; ++i)
617 auto &
min = ZoneBox[i];
618 auto &
max = ZoneBox[i + 3];
621 ZoneLength[i] =
max -
min;
624 zone->UpdateZoneBox(ZoneBox, ZoneLength);
#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.
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.
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.
std::map< int, std::map< int, CompositeSharedPtr > > & GetDomain()
void GetCompositeList(const std::string &compositeStr, CompositeMap &compositeVector) const
int GetSpaceDimension()
Dimension of the space (can be a 1D curve in 3D space).
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
std::map< int, ZoneBaseShPtr > m_zones
Movement()=default
Default constructor.
const char *const ShapeTypeMap[SIZE_ShapeType]
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
static std::string StripParentheses(const std::string &str)
std::shared_ptr< ZoneBase > ZoneBaseShPtr
std::shared_ptr< ZoneTranslate > ZoneTranslateShPtr
std::shared_ptr< InterfacePair > InterfacePairShPtr
std::shared_ptr< ZoneRotate > ZoneRotateShPtr
std::shared_ptr< Interface > InterfaceShPtr
MovementType
Enum of zone movement type.
std::shared_ptr< ZoneFixed > ZoneFixedShPtr
static std::string ReadTag(std::string &tagStr)
const std::string MovementTypeStr[]
Map of zone movement type to movement type string.
std::map< int, CompositeSharedPtr > CompositeMap
scalarT< T > max(scalarT< T > lhs, scalarT< T > rhs)
scalarT< T > min(scalarT< T > lhs, scalarT< T > rhs)