44namespace SpatialDomains
47std::string
static inline ReadTag(std::string &tagStr)
49 std::string::size_type indxBeg = tagStr.find_first_of(
'[') + 1;
50 std::string::size_type indxEnd = tagStr.find_last_of(
']') - 1;
54 (std::string(
"Error reading interface region definition:") + tagStr)
57 std::string indxStr = tagStr.substr(indxBeg, indxEnd - indxBeg + 1);
65 TiXmlNode *nektar = pSession->GetElement(
"NEKTAR");
66 if (nektar ==
nullptr)
71 TiXmlNode *movement = nektar->FirstChild(
"MOVEMENT");
72 if (movement !=
nullptr)
74 bool zones = movement->FirstChild(
"ZONES") !=
nullptr;
77 TiXmlElement *zonesTag =
78 pSession->GetElement(
"NEKTAR/MOVEMENT/ZONES");
82 bool interfaces = movement->FirstChild(
"INTERFACES") !=
nullptr;
85 TiXmlElement *interfacesTag =
86 pSession->GetElement(
"NEKTAR/MOVEMENT/INTERFACES");
91 "Only one of ZONES or INTERFACES present in the MOVEMENT "
95 if (pSession->DefinesSolverInfo(
"DiffusionType"))
97 ASSERTL0(pSession->GetSolverInfo(
"DiffusionType") ==
"LDGNS",
98 "Only LDGNS is supported as the DiffusionType in "
99 "SOLVERINFO when a MOVEMENT block is defined.")
104 if (movement !=
nullptr && pSession->DefinesCmdLineArgument(
"verbose"))
106 auto comm = pSession->GetComm();
107 if (comm->TreatAsRankZero())
109 std::cout <<
"Movement Info:\n";
110 std::cout <<
"\tNum zones: " <<
m_zones.size() <<
"\n";
115 auto numEl = zone.second->GetElements().size();
120 !zone.second->GetElements().empty()
121 ? zone.second->GetElements().front()->GetShapeType()
125 if (comm->TreatAsRankZero())
127 std::cout <<
"\t- " << zone.first <<
" "
129 zone.second->GetMovementType())]
130 <<
": " << numEl <<
" "
135 if (comm->TreatAsRankZero())
137 std::cout <<
"\tNum interfaces: " <<
m_interfaces.size() <<
"\n";
143 interface.second->GetLeftInterface()->GetEdge().size();
145 interface.second->GetRightInterface()->GetEdge().size();
151 !interface.second->GetLeftInterface()->GetEdge().empty()
152 ? interface.second->GetLeftInterface()
155 ->second->GetShapeType()
157 comm->GetSpaceComm()->AllReduce(shapeTypeLeft,
160 !interface.second->GetRightInterface()->GetEdge().empty()
161 ? interface.second->GetRightInterface()
164 ->second->GetShapeType()
166 comm->GetSpaceComm()->AllReduce(shapeTypeRight,
169 if (comm->TreatAsRankZero())
171 std::cout <<
"\t- \"" << interface.first.second <<
"\": "
172 << interface.second->GetLeftInterface()->GetId()
173 <<
" (" << numLeft <<
" "
176 << interface.second->GetRightInterface()->GetId()
177 <<
" (" << numRight <<
" "
183 comm->GetSpaceComm()->Block();
184 if (comm->TreatAsRankZero())
186 std::cout << std::endl;
196 ASSERTL0(zonesTag,
"Unable to find ZONES tag in file.");
197 TiXmlElement *zonesElement = zonesTag->FirstChildElement();
200 std::string zoneType = zonesElement->Value();
205 err = zonesElement->QueryIntAttribute(
"ID", &indx);
206 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read zone ID.");
208 std::string interfaceDomainStr;
209 err = zonesElement->QueryStringAttribute(
"DOMAIN", &interfaceDomainStr);
210 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read zone domain.");
213 auto domFind = stoi(
ReadTag(interfaceDomainStr));
214 std::map<int, CompositeSharedPtr> domain;
215 if (domains.find(domFind) != domains.end())
217 domain = domains.at(domFind);
221 if (zoneType ==
"F" || zoneType ==
"FIXED")
224 indx, domain, coordDim));
226 else if (zoneType ==
"R" || zoneType ==
"ROTATE" ||
227 zoneType ==
"ROTATING")
229 std::string originStr;
230 err = zonesElement->QueryStringAttribute(
"ORIGIN", &originStr);
231 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read origin.");
232 std::vector<NekDouble> originVec;
238 err = zonesElement->QueryStringAttribute(
"AXIS", &axisStr);
239 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read axis.");
243 std::string angularVelStr;
244 err = zonesElement->QueryStringAttribute(
"ANGVEL", &angularVelStr);
245 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read angular velocity.");
249 pSession->GetInterpreter(), angularVelStr);
252 indx, domain, coordDim, origin, axis, angularVelEqn));
256 else if (zoneType ==
"T" || zoneType ==
"TRANSLATE" ||
257 zoneType ==
"TRANSLATING")
259 std::string velocityStr;
260 err = zonesElement->QueryStringAttribute(
"VELOCITY", &velocityStr);
261 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read direction.");
262 std::vector<NekDouble> velocity;
267 indx, domain, coordDim, velocity));
271 else if (zoneType ==
"P" || zoneType ==
"PRESCRIBED")
273 std::string xDeformStr;
274 err = zonesElement->QueryStringAttribute(
"XDEFORM", &xDeformStr);
275 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read x deform equation.");
277 std::make_shared<LibUtilities::Equation>(
278 pSession->GetInterpreter(), xDeformStr);
280 std::string yDeformStr;
281 err = zonesElement->QueryStringAttribute(
"YDEFORM", &yDeformStr);
282 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read y deform equation.");
284 std::make_shared<LibUtilities::Equation>(
285 pSession->GetInterpreter(), yDeformStr);
287 std::string zDeformStr;
288 err = zonesElement->QueryStringAttribute(
"ZDEFORM", &zDeformStr);
289 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read z deform equation.");
291 std::make_shared<LibUtilities::Equation>(
292 pSession->GetInterpreter(), zDeformStr);
296 indx, domain, coordDim, xDeformEqn, yDeformEqn,
303 WARNINGL0(
false,
"Zone type '" + zoneType +
304 "' is unsupported. Valid types are: 'Fixed', "
305 "'Rotate', 'Translate', or 'Prescribe'.")
309 zonesElement = zonesElement->NextSiblingElement();
315 ASSERTL0(interfacesTag,
"Unable to find INTERFACES tag in file.");
316 TiXmlElement *interfaceElement = interfacesTag->FirstChildElement();
318 while (interfaceElement)
321 "INTERFACE" == (std::string)interfaceElement->Value(),
322 "Only INTERFACE tags may be present inside the INTERFACES block.")
327 err = interfaceElement->QueryStringAttribute(
"NAME", &
name);
328 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read interface name.");
329 TiXmlElement *sideElement = interfaceElement->FirstChildElement();
335 std::vector<int> cnt;
340 "In INTERFACE NAME " +
name +
341 ", only two sides may be present in each INTERFACE block.")
344 err = sideElement->QueryIntAttribute(
"ID", &indx);
345 ASSERTL0(err == TIXML_SUCCESS,
"Unable to read interface ID.");
347 std::string boundaryStr;
349 sideElement->QueryStringAttribute(
"BOUNDARY", &boundaryStr);
352 if (boundaryErr == TIXML_SUCCESS)
354 std::string indxStr =
ReadTag(boundaryStr);
359 auto sideElVal = sideElement->ValueStr();
360 if (sideElVal ==
"LEFT" || sideElVal ==
"L")
364 else if (sideElVal ==
"RIGHT" || sideElVal ==
"R")
371 sideElement->ValueStr() +
372 " is not a valid interface side for interface "
374 name +
". Please only use LEFT or RIGHT.")
377 interfaces[cnt[cnt.size() - 1]] =
379 indx, boundaryEdge));
381 sideElement = sideElement->NextSiblingElement();
384 ASSERTL0(std::accumulate(cnt.begin(), cnt.end(), 0) == 1,
385 "You must have only one LEFT and one RIGHT side"
386 " present in interface NAME " +
391 interfaces[0], interfaces[1]));
392 interfaceElement = interfaceElement->NextSiblingElement();
400 std::set<int> movedZoneIds;
403 if (zone.second->Move(timeStep))
405 movedZoneIds.insert(zone.first);
412 int leftId = interPair.second->GetLeftInterface()->GetId();
413 int rightId = interPair.second->GetRightInterface()->GetId();
415 if (movedZoneIds.find(leftId) != movedZoneIds.end() ||
416 movedZoneIds.find(rightId) != movedZoneIds.end())
418 m_zones[leftId]->GetMoved() =
true;
419 m_zones[rightId]->GetMoved() =
true;
#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)
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 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 PerformMovement(NekDouble timeStep)
Movement(const LibUtilities::SessionReaderSharedPtr &pSession, MeshGraph *meshGraph)
Constructor.
void ReadZones(TiXmlElement *zonesTag, MeshGraph *meshGraph, const LibUtilities::SessionReaderSharedPtr &pSession)
Read zones given TiXmlDocument.
void ReadInterfaces(TiXmlElement *interfacesTag, MeshGraph *meshGraph)
Read interfaces given TiXmlDocument.
InterfaceCollection m_interfaces
std::map< int, ZoneBaseShPtr > m_zones
const char *const ShapeTypeMap[SIZE_ShapeType]
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
std::shared_ptr< ZoneBase > ZoneBaseShPtr
std::shared_ptr< ZonePrescribe > ZonePrescribeShPtr
std::shared_ptr< ZoneTranslate > ZoneTranslateShPtr
std::shared_ptr< InterfacePair > InterfacePairShPtr
std::shared_ptr< ZoneRotate > ZoneRotateShPtr
std::shared_ptr< Interface > InterfaceShPtr
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
The above copyright notice and this permission notice shall be included.