Nektar++
Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | Private Member Functions | List of all members
Nektar::SpatialDomains::Movement Class Reference

#include <Movement.h>

Public Member Functions

 Movement ()=default
 Default constructor.
 
 Movement (const LibUtilities::SessionReaderSharedPtr &pSession, MeshGraph *meshGraph)
 Constructor to read from XML file.
 
 ~Movement ()=default
 Default destructor.
 
void WriteMovement (TiXmlElement *root)
 Write the MOVEMENT section of the XML file.
 
const InterfaceCollectionGetInterfaces () const
 
const std::map< int, ZoneBaseShPtr > & GetZones () const
 
void PerformMovement (NekDouble timeStep)
 
const bool & GetMoveFlag () const
 
bool & GetCoordExchangeFlag ()
 
const bool & GetMovedFlag () const
 
const bool & GetTranslateFlag () const
 
const bool & GetRotateFlag () const
 
const bool & GetSectorRotateFlag () const
 
const bool & GetImplicitALESolverFlag () const
 
const bool & GetMeshDistortedFlag () const
 
void SetImplicitALEFlag (bool &ImplicitALE)
 
void SetMeshDistortedFlag (bool &meshDistorted)
 
void AddZone (ZoneBaseShPtr zone)
 Add a zone object to this Movement data.
 
void AddInterface (std::string name, InterfaceShPtr left, InterfaceShPtr right)
 Add pair of interfaces to this data.
 

Protected Attributes

InterfaceCollection m_interfaces
 
std::map< int, ZoneBaseShPtrm_zones
 
bool m_moveFlag = false
 
bool m_translate = false
 
bool m_rotate = false
 
bool m_sectorRotate = false
 
bool m_moved = false
 
bool m_implicitALESolver = false
 
bool m_meshDistorted = false
 
bool m_coordExchangeFlag
 

Private Member Functions

void ReadZones (TiXmlElement *zonesTag, MeshGraph *meshGraph, const LibUtilities::SessionReaderSharedPtr &pSession)
 Read zones given TiXmlDocument.
 
void ReadInterfaces (TiXmlElement *interfacesTag, MeshGraph *meshGraph)
 Read interfaces given TiXmlDocument.
 
void UpdateTransZoneBox (const LibUtilities::SessionReaderSharedPtr &pSession)
 Update the box of translate zones.
 

Detailed Description

Definition at line 52 of file Movement.h.

Constructor & Destructor Documentation

◆ Movement() [1/2]

Nektar::SpatialDomains::Movement::Movement ( )
default

Default constructor.

◆ Movement() [2/2]

Nektar::SpatialDomains::Movement::Movement ( const LibUtilities::SessionReaderSharedPtr pSession,
MeshGraph meshGraph 
)

Constructor to read from XML file.

Definition at line 71 of file Movement/Movement.cpp.

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}
void UpdateTransZoneBox(const LibUtilities::SessionReaderSharedPtr &pSession)
Update the box of translate zones.
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
Definition Movement.h:141
std::map< int, ZoneBaseShPtr > m_zones
Definition Movement.h:142
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition ShapeType.hpp:75
const std::string MovementTypeStr[]
Map of zone movement type to movement type string.
Definition Zones.h:57

References m_interfaces, m_translate, m_zones, Nektar::SpatialDomains::MovementTypeStr, ReadInterfaces(), ReadZones(), Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::ReduceSum, Nektar::LibUtilities::ShapeTypeMap, and UpdateTransZoneBox().

◆ ~Movement()

Nektar::SpatialDomains::Movement::~Movement ( )
default

Default destructor.

Member Function Documentation

◆ AddInterface()

void Nektar::SpatialDomains::Movement::AddInterface ( std::string  name,
InterfaceShPtr  left,
InterfaceShPtr  right 
)

Add pair of interfaces to this data.

Store an interface pair with this Movement data.

Definition at line 592 of file Movement/Movement.cpp.

594{
595 m_interfaces[std::make_pair(m_interfaces.size(), name)] =
598}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< InterfacePair > InterfacePairShPtr

References m_interfaces.

Referenced by Nektar::MovementTests::BOOST_AUTO_TEST_CASE(), Nektar::MovementTests::BOOST_AUTO_TEST_CASE(), and export_Movement().

◆ AddZone()

void Nektar::SpatialDomains::Movement::AddZone ( ZoneBaseShPtr  zone)

Add a zone object to this Movement data.

Store a zone object with this Movement data.

Definition at line 581 of file Movement/Movement.cpp.

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}
MovementType
Enum of zone movement type.
Definition Zones.h:48

References Nektar::SpatialDomains::eFixed, Nektar::SpatialDomains::eNone, m_moveFlag, and m_zones.

Referenced by Nektar::MovementTests::BOOST_AUTO_TEST_CASE(), Nektar::MovementTests::BOOST_AUTO_TEST_CASE(), and export_Movement().

◆ GetCoordExchangeFlag()

bool & Nektar::SpatialDomains::Movement::GetCoordExchangeFlag ( )
inline

Definition at line 86 of file Movement.h.

87 {
89 }

References m_coordExchangeFlag.

◆ GetImplicitALESolverFlag()

const bool & Nektar::SpatialDomains::Movement::GetImplicitALESolverFlag ( ) const
inline

Definition at line 111 of file Movement.h.

112 {
113 return m_implicitALESolver;
114 }

References m_implicitALESolver.

◆ GetInterfaces()

const InterfaceCollection & Nektar::SpatialDomains::Movement::GetInterfaces ( ) const
inline

Definition at line 69 of file Movement.h.

70 {
71 return m_interfaces;
72 }

References m_interfaces.

Referenced by Nektar::MovementTests::BOOST_AUTO_TEST_CASE(), and Nektar::MovementTests::BOOST_AUTO_TEST_CASE().

◆ GetMeshDistortedFlag()

const bool & Nektar::SpatialDomains::Movement::GetMeshDistortedFlag ( ) const
inline

Definition at line 116 of file Movement.h.

117 {
118 return m_meshDistorted;
119 }

References m_meshDistorted.

◆ GetMovedFlag()

const bool & Nektar::SpatialDomains::Movement::GetMovedFlag ( ) const
inline

Definition at line 91 of file Movement.h.

92 {
93 return m_moved;
94 }

References m_moved.

◆ GetMoveFlag()

const bool & Nektar::SpatialDomains::Movement::GetMoveFlag ( ) const
inline

Definition at line 81 of file Movement.h.

82 {
83 return m_moveFlag;
84 }

References m_moveFlag.

◆ GetRotateFlag()

const bool & Nektar::SpatialDomains::Movement::GetRotateFlag ( ) const
inline

Definition at line 101 of file Movement.h.

102 {
103 return m_rotate;
104 }

References m_rotate.

◆ GetSectorRotateFlag()

const bool & Nektar::SpatialDomains::Movement::GetSectorRotateFlag ( ) const
inline

Definition at line 106 of file Movement.h.

107 {
108 return m_sectorRotate;
109 }

References m_sectorRotate.

◆ GetTranslateFlag()

const bool & Nektar::SpatialDomains::Movement::GetTranslateFlag ( ) const
inline

Definition at line 96 of file Movement.h.

97 {
98 return m_translate;
99 }

References m_translate.

◆ GetZones()

const std::map< int, ZoneBaseShPtr > & Nektar::SpatialDomains::Movement::GetZones ( ) const
inline

Definition at line 74 of file Movement.h.

75 {
76 return m_zones;
77 }

References m_zones.

Referenced by Nektar::MovementTests::BOOST_AUTO_TEST_CASE(), Nektar::MovementTests::BOOST_AUTO_TEST_CASE(), and export_Movement().

◆ PerformMovement()

void Nektar::SpatialDomains::Movement::PerformMovement ( NekDouble  timeStep)

Definition at line 553 of file Movement/Movement.cpp.

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}

References m_interfaces, m_moved, and m_zones.

Referenced by export_Movement().

◆ ReadInterfaces()

void Nektar::SpatialDomains::Movement::ReadInterfaces ( TiXmlElement *  interfacesTag,
MeshGraph meshGraph 
)
private

Read interfaces given TiXmlDocument.

Definition at line 368 of file Movement/Movement.cpp.

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}
NekDouble L
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
std::shared_ptr< Interface > InterfaceShPtr
static std::string ReadTag(std::string &tagStr)
std::map< int, CompositeSharedPtr > CompositeMap
Definition MeshGraph.h:179
STL namespace.

References ASSERTL0, Nektar::ErrorUtil::efatal, Nektar::SpatialDomains::MeshGraph::GetCompositeList(), m_interfaces, NEKERROR, and Nektar::SpatialDomains::ReadTag().

Referenced by Movement().

◆ ReadZones()

void Nektar::SpatialDomains::Movement::ReadZones ( TiXmlElement *  zonesTag,
MeshGraph meshGraph,
const LibUtilities::SessionReaderSharedPtr pSession 
)
private

Read zones given TiXmlDocument.

Definition at line 193 of file Movement/Movement.cpp.

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;
266 Array<OneD, NekDouble> base(3);
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
314 Array<OneD, LibUtilities::EquationSharedPtr> velocityEqns(coordDim);
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
339 Array<OneD, LibUtilities::EquationSharedPtr> displacementEqns(
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}
#define WARNINGL0(condition, msg)
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
std::shared_ptr< Equation > EquationSharedPtr
Definition Equation.h:131
const NekDouble sector
std::vector< std::string > velocityStr
const NekPoint< NekDouble > origin
const NekDouble rampTime
const Array< OneD, NekDouble > base(3, 0.0)
std::vector< std::string > displacementStr
std::shared_ptr< ZoneBase > ZoneBaseShPtr
Definition Zones.h:177
std::shared_ptr< ZoneTranslate > ZoneTranslateShPtr
Definition Zones.h:376
std::shared_ptr< ZoneRotate > ZoneRotateShPtr
Definition Zones.h:375
std::shared_ptr< ZoneFixed > ZoneFixedShPtr
Definition Zones.h:377
NekVector< NekDouble > DNekVec

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::Interpreter::DefineFunction(), Nektar::LibUtilities::Interpreter::Evaluate(), Nektar::ParseUtils::GenerateVector(), Nektar::SpatialDomains::MeshGraph::GetDomain(), Nektar::SpatialDomains::MeshGraph::GetSpaceDimension(), m_moveFlag, m_rotate, m_sectorRotate, m_translate, m_zones, Nektar::NekVector< DataType >::Normalize(), Nektar::SpatialDomains::ReadTag(), and WARNINGL0.

Referenced by Movement().

◆ SetImplicitALEFlag()

void Nektar::SpatialDomains::Movement::SetImplicitALEFlag ( bool &  ImplicitALE)
inline

Definition at line 121 of file Movement.h.

122 {
123 m_implicitALESolver = ImplicitALE;
124 }

References m_implicitALESolver.

◆ SetMeshDistortedFlag()

void Nektar::SpatialDomains::Movement::SetMeshDistortedFlag ( bool &  meshDistorted)
inline

Definition at line 126 of file Movement.h.

127 {
128 m_moved = meshDistorted;
129 }

References m_moved.

◆ UpdateTransZoneBox()

void Nektar::SpatialDomains::Movement::UpdateTransZoneBox ( const LibUtilities::SessionReaderSharedPtr pSession)
private

Update the box of translate zones.

Generate domain box for translation mesh.

Definition at line 601 of file Movement/Movement.cpp.

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}
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

References Nektar::SpatialDomains::eTranslate, m_zones, tinysimd::max(), tinysimd::min(), Nektar::LibUtilities::ReduceMax, and Nektar::LibUtilities::ReduceMin.

Referenced by Movement().

◆ WriteMovement()

void Nektar::SpatialDomains::Movement::WriteMovement ( TiXmlElement *  root)

Write the MOVEMENT section of the XML file.

Export this Movement information to a Nektar++ XML file.

Definition at line 454 of file Movement/Movement.cpp.

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}
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 std::string StripParentheses(const std::string &str)
std::vector< double > z(NPUPPER)

References Nektar::SpatialDomains::eRotate, Nektar::SpatialDomains::eTranslate, Nektar::ParseUtils::GenerateSeqString(), m_interfaces, m_zones, Nektar::SpatialDomains::MovementTypeStr, and Nektar::SpatialDomains::StripParentheses().

Referenced by Nektar::MovementTests::BOOST_AUTO_TEST_CASE().

Member Data Documentation

◆ m_coordExchangeFlag

bool Nektar::SpatialDomains::Movement::m_coordExchangeFlag
protected
Initial value:
=
true

Definition at line 150 of file Movement.h.

Referenced by GetCoordExchangeFlag().

◆ m_implicitALESolver

bool Nektar::SpatialDomains::Movement::m_implicitALESolver = false
protected

Definition at line 148 of file Movement.h.

Referenced by GetImplicitALESolverFlag(), and SetImplicitALEFlag().

◆ m_interfaces

InterfaceCollection Nektar::SpatialDomains::Movement::m_interfaces
protected

◆ m_meshDistorted

bool Nektar::SpatialDomains::Movement::m_meshDistorted = false
protected

Definition at line 149 of file Movement.h.

Referenced by GetMeshDistortedFlag().

◆ m_moved

bool Nektar::SpatialDomains::Movement::m_moved = false
protected

Definition at line 147 of file Movement.h.

Referenced by GetMovedFlag(), PerformMovement(), and SetMeshDistortedFlag().

◆ m_moveFlag

bool Nektar::SpatialDomains::Movement::m_moveFlag = false
protected

Definition at line 143 of file Movement.h.

Referenced by AddZone(), GetMoveFlag(), and ReadZones().

◆ m_rotate

bool Nektar::SpatialDomains::Movement::m_rotate = false
protected

Definition at line 145 of file Movement.h.

Referenced by GetRotateFlag(), and ReadZones().

◆ m_sectorRotate

bool Nektar::SpatialDomains::Movement::m_sectorRotate = false
protected

Definition at line 146 of file Movement.h.

Referenced by GetSectorRotateFlag(), and ReadZones().

◆ m_translate

bool Nektar::SpatialDomains::Movement::m_translate = false
protected

Definition at line 144 of file Movement.h.

Referenced by GetTranslateFlag(), Movement(), and ReadZones().

◆ m_zones

std::map<int, ZoneBaseShPtr> Nektar::SpatialDomains::Movement::m_zones
protected