Nektar++
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. More...
 
 Movement (const LibUtilities::SessionReaderSharedPtr &pSession, MeshGraph *meshGraph)
 Constructor to read from XML file. More...
 
 ~Movement ()=default
 Default destructor. More...
 
void WriteMovement (TiXmlElement *root)
 Write the MOVEMENT section of the XML file. More...
 
const InterfaceCollectionGetInterfaces () const
 
const std::map< int, ZoneBaseShPtr > & GetZones () const
 
void PerformMovement (NekDouble timeStep)
 
const bool & GetMoveFlag () const
 
bool & GetCoordExchangeFlag ()
 
const Array< OneD, NekDouble > & GetDomainBox () const
 
const Array< OneD, NekDouble > & GetDomainLength () const
 
const bool & GetMovedFlag () const
 
const bool & GetTranslateFlag () const
 
const bool & GetImplicitALESolverFlag () const
 
void SetImplicitALEFlag (bool &ImplicitALE)
 
void AddZone (ZoneBaseShPtr zone)
 Add a zone object to this Movement data. More...
 
void AddInterface (std::string name, InterfaceShPtr left, InterfaceShPtr right)
 Add pair of interfaces to this data. More...
 

Protected Attributes

InterfaceCollection m_interfaces
 
std::map< int, ZoneBaseShPtrm_zones
 
bool m_moveFlag = false
 
bool m_translate = false
 
bool m_moved = false
 
bool m_ImplicitALESolver = false
 
bool m_coordExchangeFlag
 
Array< OneD, NekDoublem_DomainBox
 
Array< OneD, NekDoublem_DomainLength
 

Private Member Functions

void ReadZones (TiXmlElement *zonesTag, MeshGraph *meshGraph, const LibUtilities::SessionReaderSharedPtr &pSession)
 Read zones given TiXmlDocument. More...
 
void ReadInterfaces (TiXmlElement *interfacesTag, MeshGraph *meshGraph)
 Read interfaces given TiXmlDocument. More...
 
void DomainBox ()
 Calculate length of the domain. More...
 

Detailed Description

Definition at line 51 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 70 of file Movement/Movement.cpp.

72{
73 TiXmlNode *nektar = pSession->GetElement("NEKTAR");
74 if (nektar == nullptr)
75 {
76 return;
77 }
78
79 TiXmlNode *movement = nektar->FirstChild("MOVEMENT");
80 if (movement != nullptr)
81 {
82 bool zones = movement->FirstChild("ZONES") != nullptr;
83 if (zones)
84 {
85 TiXmlElement *zonesTag =
86 pSession->GetElement("NEKTAR/MOVEMENT/ZONES");
87 ReadZones(zonesTag, meshGraph, pSession);
88 }
89
90 bool interfaces = movement->FirstChild("INTERFACES") != nullptr;
91 if (interfaces)
92 {
93 TiXmlElement *interfacesTag =
94 pSession->GetElement("NEKTAR/MOVEMENT/INTERFACES");
95 ReadInterfaces(interfacesTag, meshGraph);
96 }
97
98 ASSERTL0(zones == interfaces,
99 "Only one of ZONES or INTERFACES present in the MOVEMENT "
100 "block.")
101 // Generate domain box
102 DomainBox();
103 // Get DomainBox from all processes
104 if (m_translate)
105 {
106 auto comm = pSession->GetComm();
107 for (int i = 0; i < 3; ++i)
108 {
109 auto &min = m_DomainBox[i];
110 auto &max = m_DomainBox[i + 3];
111 comm->GetSpaceComm()->AllReduce(min, LibUtilities::ReduceMin);
112 comm->GetSpaceComm()->AllReduce(max, LibUtilities::ReduceMax);
113 m_DomainLength[i] = max - min;
114 }
115 }
116 }
117
118 // DEBUG COMMENTS
119 if (movement != nullptr && pSession->DefinesCmdLineArgument("verbose"))
120 {
121 auto comm = pSession->GetComm();
122 if (comm->TreatAsRankZero())
123 {
124 std::cout << "Movement Info:\n";
125 std::cout << "\tNum zones: " << m_zones.size() << "\n";
126 }
127
128 for (auto &zone : m_zones)
129 {
130 auto numEl = zone.second->GetElements().size();
131 comm->GetSpaceComm()->AllReduce(numEl, LibUtilities::ReduceSum);
132
133 // Find shape type if not on this proc
134 int shapeType =
135 !zone.second->GetElements().empty()
136 ? zone.second->GetElements().front()->GetShapeType()
137 : -1;
138 comm->GetSpaceComm()->AllReduce(shapeType, LibUtilities::ReduceMax);
139
140 if (comm->TreatAsRankZero())
141 {
142 std::cout << "\t- " << zone.first << " "
143 << MovementTypeStr[static_cast<int>(
144 zone.second->GetMovementType())]
145 << ": " << numEl << " "
146 << LibUtilities::ShapeTypeMap[shapeType] << "s\n";
147 }
148 }
149
150 if (comm->TreatAsRankZero())
151 {
152 std::cout << "\tNum interfaces: " << m_interfaces.size() << "\n";
153 }
154
155 for (auto &interface : m_interfaces)
156 {
157 auto numLeft =
158 interface.second->GetLeftInterface()->GetEdge().size();
159 auto numRight =
160 interface.second->GetRightInterface()->GetEdge().size();
161 comm->GetSpaceComm()->AllReduce(numLeft, LibUtilities::ReduceSum);
162 comm->GetSpaceComm()->AllReduce(numRight, LibUtilities::ReduceSum);
163
164 // Find shape type if not on this proc
165 int shapeTypeLeft =
166 !interface.second->GetLeftInterface()->GetEdge().empty()
167 ? interface.second->GetLeftInterface()
168 ->GetEdge()
169 .begin()
170 ->second->GetShapeType()
171 : -1;
172 comm->GetSpaceComm()->AllReduce(shapeTypeLeft,
174 int shapeTypeRight =
175 !interface.second->GetRightInterface()->GetEdge().empty()
176 ? interface.second->GetRightInterface()
177 ->GetEdge()
178 .begin()
179 ->second->GetShapeType()
180 : -1;
181 comm->GetSpaceComm()->AllReduce(shapeTypeRight,
183
184 if (comm->TreatAsRankZero())
185 {
186 std::cout << "\t- \"" << interface.first.second << "\": "
187 << interface.second->GetLeftInterface()->GetId()
188 << " (" << numLeft << " "
189 << LibUtilities::ShapeTypeMap[shapeTypeLeft]
190 << "s) <-> "
191 << interface.second->GetRightInterface()->GetId()
192 << " (" << numRight << " "
193 << LibUtilities::ShapeTypeMap[shapeTypeRight]
194 << "s)\n";
195 }
196 }
197
198 comm->GetSpaceComm()->Block();
199 if (comm->TreatAsRankZero())
200 {
201 std::cout << std::endl;
202 }
203 }
204}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void DomainBox()
Calculate length of the domain.
Array< OneD, NekDouble > m_DomainLength
Definition: Movement.h:139
void ReadZones(TiXmlElement *zonesTag, MeshGraph *meshGraph, const LibUtilities::SessionReaderSharedPtr &pSession)
Read zones given TiXmlDocument.
void ReadInterfaces(TiXmlElement *interfacesTag, MeshGraph *meshGraph)
Read interfaces given TiXmlDocument.
Array< OneD, NekDouble > m_DomainBox
Definition: Movement.h:138
InterfaceCollection m_interfaces
Definition: Movement.h:130
std::map< int, ZoneBaseShPtr > m_zones
Definition: Movement.h:131
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:56

References ASSERTL0, DomainBox(), m_DomainBox, m_DomainLength, m_interfaces, m_translate, m_zones, Nektar::SpatialDomains::MovementTypeStr, ReadInterfaces(), ReadZones(), Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::ReduceMin, Nektar::LibUtilities::ReduceSum, and Nektar::LibUtilities::ShapeTypeMap.

◆ ~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 567 of file Movement/Movement.cpp.

569{
570 m_interfaces[std::make_pair(m_interfaces.size(), name)] =
573}
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, and CellMLToNektar.pycml::name.

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

◆ 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 556 of file Movement/Movement.cpp.

557{
558 m_zones[zone->GetId()] = zone;
559 MovementType mtype = zone->GetMovementType();
560 if (mtype != MovementType::eFixed && mtype != MovementType::eNone)
561 {
562 m_moveFlag = true;
563 }
564}
MovementType
Enum of zone movement type.
Definition: Zones.h:47

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

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

◆ DomainBox()

void Nektar::SpatialDomains::Movement::DomainBox ( )
private

Calculate length of the domain.

Generate domain box for translation mesh.

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

577{
578 if (m_translate)
579 {
580 // NekDouble minx, miny, minz, maxx, maxy, maxz;
581 Array<OneD, NekDouble> min(3, 0.0), max(3, 0.0);
582 Array<OneD, NekDouble> x(3, 0.0);
583
584 // Find a point in the domian for parallel
585 for (auto &zones : m_zones)
586 {
587 int NumVerts = zones.second->GetOriginalVertex().size();
588 if (NumVerts != 0)
589 {
590 PointGeom p = zones.second->GetOriginalVertex()[0];
591 p.GetCoords(x[0], x[1], x[2]);
592 for (int j = 0; j < 3; ++j)
593 {
594 min[j] = x[j];
595 max[j] = x[j];
596 }
597 break;
598 }
599 }
600
601 // loop over all zones and original vertexes
602 for (auto &zones : m_zones)
603 {
604 int NumVerts = zones.second->GetOriginalVertex().size();
605 for (int i = 0; i < NumVerts; ++i)
606 {
607 PointGeom p = zones.second->GetOriginalVertex()[i];
608 p.GetCoords(x[0], x[1], x[2]);
609 for (int j = 0; j < 3; ++j)
610 {
611 min[j] = (x[j] < min[j] ? x[j] : min[j]);
612 max[j] = (x[j] > max[j] ? x[j] : max[j]);
613 }
614 }
615 }
616 // save bounding length
617 m_DomainLength = Array<OneD, NekDouble>(3, 0.0);
618 // save bounding box
619 m_DomainBox = Array<OneD, NekDouble>(6, 0.0);
620 for (int j = 0; j < 3; ++j)
621 {
622 m_DomainBox[j] = min[j];
623 m_DomainBox[j + 3] = max[j];
624 }
625 }
626}

References m_DomainBox, m_DomainLength, m_translate, m_zones, and CellMLToNektar.cellml_metadata::p.

Referenced by Movement().

◆ GetCoordExchangeFlag()

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

Definition at line 85 of file Movement.h.

86 {
88 }

References m_coordExchangeFlag.

◆ GetDomainBox()

const Array< OneD, NekDouble > & Nektar::SpatialDomains::Movement::GetDomainBox ( ) const
inline

Definition at line 90 of file Movement.h.

91 {
92 return m_DomainBox;
93 }

References m_DomainBox.

◆ GetDomainLength()

const Array< OneD, NekDouble > & Nektar::SpatialDomains::Movement::GetDomainLength ( ) const
inline

Definition at line 95 of file Movement.h.

96 {
97 return m_DomainLength;
98 }

References m_DomainLength.

◆ GetImplicitALESolverFlag()

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

Definition at line 110 of file Movement.h.

111 {
112 return m_ImplicitALESolver;
113 }

References m_ImplicitALESolver.

◆ GetInterfaces()

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

Definition at line 68 of file Movement.h.

69 {
70 return m_interfaces;
71 }

References m_interfaces.

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

◆ GetMovedFlag()

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

Definition at line 100 of file Movement.h.

101 {
102 return m_moved;
103 }

References m_moved.

◆ GetMoveFlag()

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

Definition at line 80 of file Movement.h.

81 {
82 return m_moveFlag;
83 }

References m_moveFlag.

◆ GetTranslateFlag()

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

Definition at line 105 of file Movement.h.

106 {
107 return m_translate;
108 }

References m_translate.

◆ GetZones()

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

Definition at line 73 of file Movement.h.

74 {
75 return m_zones;
76 }

References m_zones.

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

◆ PerformMovement()

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

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

529{
530 std::set<int> movedZoneIds;
531 for (auto &zone : m_zones)
532 {
533 if (zone.second->Move(timeStep))
534 {
535 movedZoneIds.insert(zone.first);
536 }
537 }
538
539 // If zone has moved, set all interfaces on that zone to moved.
540 for (auto &interPair : m_interfaces)
541 {
542 int leftId = interPair.second->GetLeftInterface()->GetId();
543 int rightId = interPair.second->GetRightInterface()->GetId();
544
545 if (movedZoneIds.find(leftId) != movedZoneIds.end() ||
546 movedZoneIds.find(rightId) != movedZoneIds.end())
547 {
548 m_zones[leftId]->GetMoved() = true;
549 m_zones[rightId]->GetMoved() = true;
550 }
551 }
552 m_moved = true;
553}

References m_interfaces, m_moved, and m_zones.

◆ ReadInterfaces()

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

Read interfaces given TiXmlDocument.

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

344{
345 ASSERTL0(interfacesTag, "Unable to find INTERFACES tag in file.");
346 TiXmlElement *interfaceElement = interfacesTag->FirstChildElement();
347
348 while (interfaceElement)
349 {
350 ASSERTL0(
351 "INTERFACE" == (std::string)interfaceElement->Value(),
352 "Only INTERFACE tags may be present inside the INTERFACES block.")
353
354 int err;
355
356 std::string name;
357 err = interfaceElement->QueryStringAttribute("NAME", &name);
358 ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface name.");
359 TiXmlElement *sideElement = interfaceElement->FirstChildElement();
360
361 Array<OneD, InterfaceShPtr> interfaces(2);
362 std::vector<int> cnt;
363 while (sideElement)
364 {
365 ASSERTL0(
366 cnt.size() < 2,
367 "In INTERFACE NAME " + name +
368 ", only two sides may be present in each INTERFACE block.")
369
370 // Read ID
371 int indx;
372 err = sideElement->QueryIntAttribute("ID", &indx);
373 ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface ID.");
374
375 // Read boundary
376 std::string boundaryStr;
377 err = sideElement->QueryStringAttribute("BOUNDARY", &boundaryStr);
378 ASSERTL0(err == TIXML_SUCCESS,
379 "Unable to read BOUNDARY in interface " +
380 std::to_string(indx));
381
382 CompositeMap boundaryEdge;
383 std::string indxStr = ReadTag(boundaryStr);
384 meshGraph->GetCompositeList(indxStr, boundaryEdge);
385
386 // Read SKIPCHECK flag
387 bool skipCheck = false;
388 err = sideElement->QueryBoolAttribute("SKIPCHECK", &skipCheck);
389
390 // Sets location in interface pair to 0 for left, and 1 for right
391 auto sideElVal = sideElement->ValueStr();
392 if (sideElVal == "LEFT" || sideElVal == "L")
393 {
394 cnt.emplace_back(0);
395 }
396 else if (sideElVal == "RIGHT" || sideElVal == "R")
397 {
398 cnt.emplace_back(1);
399 }
400 else
401 {
403 sideElement->ValueStr() +
404 " is not a valid interface side for interface "
405 "NAME " +
406 name + ". Please only use LEFT or RIGHT.")
407 }
408
409 interfaces[cnt[cnt.size() - 1]] =
411 indx, boundaryEdge, skipCheck));
412
413 sideElement = sideElement->NextSiblingElement();
414 }
415
416 ASSERTL0(std::accumulate(cnt.begin(), cnt.end(), 0) == 1,
417 "You must have only one LEFT and one RIGHT side"
418 " present in interface NAME " +
419 name)
420
421 m_interfaces[std::make_pair(m_interfaces.size(), name)] =
423 interfaces[0], interfaces[1]));
424 interfaceElement = interfaceElement->NextSiblingElement();
425 }
426}
NekDouble L
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
std::shared_ptr< Interface > InterfaceShPtr
static std::string ReadTag(std::string &tagStr)
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:136
STL namespace.

References ASSERTL0, Nektar::ErrorUtil::efatal, Nektar::SpatialDomains::MeshGraph::GetCompositeList(), m_interfaces, CellMLToNektar.pycml::name, 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 206 of file Movement/Movement.cpp.

208{
209 int coordDim = meshGraph->GetSpaceDimension();
210
211 ASSERTL0(zonesTag, "Unable to find ZONES tag in file.");
212 TiXmlElement *zonesElement = zonesTag->FirstChildElement();
213 while (zonesElement)
214 {
215 std::string zoneType = zonesElement->Value();
216
217 int err;
218 int indx;
219
220 err = zonesElement->QueryIntAttribute("ID", &indx);
221 ASSERTL0(err == TIXML_SUCCESS, "Unable to read zone ID.");
222
223 std::string interfaceDomainStr;
224 err = zonesElement->QueryStringAttribute("DOMAIN", &interfaceDomainStr);
225 ASSERTL0(err == TIXML_SUCCESS, "Unable to read zone domain.");
226
227 auto &domains = meshGraph->GetDomain();
228 auto domFind = std::stoi(ReadTag(interfaceDomainStr));
229 std::map<int, CompositeSharedPtr> domain;
230 if (domains.find(domFind) != domains.end())
231 {
232 domain = domains.at(domFind);
233 }
234
235 ZoneBaseShPtr zone;
236 if (zoneType == "F" || zoneType == "FIXED")
237 {
239 indx, domFind, domain, coordDim));
240 }
241 else if (zoneType == "R" || zoneType == "ROTATE" ||
242 zoneType == "ROTATING")
243 {
244 std::string originStr;
245 err = zonesElement->QueryStringAttribute("ORIGIN", &originStr);
246 ASSERTL0(err == TIXML_SUCCESS, "Unable to read origin.");
247 std::vector<NekDouble> originVec;
248 ParseUtils::GenerateVector(originStr, originVec);
249 NekPoint<NekDouble> origin =
250 NekPoint<NekDouble>(originVec[0], originVec[1], originVec[2]);
251
252 std::string axisStr;
253 err = zonesElement->QueryStringAttribute("AXIS", &axisStr);
254 ASSERTL0(err == TIXML_SUCCESS, "Unable to read axis.");
255 DNekVec axis(axisStr);
256 axis.Normalize();
257
258 std::string angularVelStr;
259 err = zonesElement->QueryStringAttribute("ANGVEL", &angularVelStr);
260 ASSERTL0(err == TIXML_SUCCESS, "Unable to read angular velocity.");
261
262 LibUtilities::EquationSharedPtr angularVelEqn =
264 pSession->GetInterpreter(), angularVelStr);
265
267 indx, domFind, domain, coordDim, origin, axis, angularVelEqn));
268
269 m_moveFlag = true;
270 }
271 else if (zoneType == "T" || zoneType == "TRANSLATE" ||
272 zoneType == "TRANSLATING")
273 {
274 // Read velocity information
275 std::string velocityStr;
276 err = zonesElement->QueryStringAttribute("VELOCITY", &velocityStr);
277 ASSERTL0(err == TIXML_SUCCESS,
278 "Unable to read VELOCITY for translating zone " +
279 std::to_string(indx));
280
281 std::vector<std::string> velocityStrSplit;
282 ParseUtils::GenerateVector(velocityStr, velocityStrSplit);
283
284 ASSERTL0(velocityStrSplit.size() >= coordDim,
285 "VELOCITY dimensions must be greater than or equal to the "
286 "coordinate dimension for translating zone " +
287 std::to_string(indx));
288
289 Array<OneD, LibUtilities::EquationSharedPtr> velocityEqns(coordDim);
290 for (int i = 0; i < coordDim; ++i)
291 {
292 velocityEqns[i] =
294 pSession->GetInterpreter(), velocityStrSplit[i]);
295 }
296
297 // Read displacement information
298 std::string displacementStr;
299 err = zonesElement->QueryStringAttribute("DISPLACEMENT",
301 ASSERTL0(err == TIXML_SUCCESS,
302 "Unable to read DISPLACEMENT for translating zone " +
303 std::to_string(indx));
304
305 std::vector<std::string> displacementStrSplit;
306 ParseUtils::GenerateVector(displacementStr, displacementStrSplit);
307
308 ASSERTL0(
309 displacementStrSplit.size() >= coordDim,
310 "DISPLACEMENT 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> displacementEqns(
315 coordDim);
316 for (int i = 0; i < coordDim; ++i)
317 {
318 displacementEqns[i] =
320 pSession->GetInterpreter(), displacementStrSplit[i]);
321 }
322
323 zone = ZoneTranslateShPtr(
325 indx, domFind, domain, coordDim, velocityEqns,
326 displacementEqns));
327
328 m_moveFlag = true;
329 m_translate = true;
330 }
331 else
332 {
333 WARNINGL0(false, "Zone type '" + zoneType +
334 "' is unsupported. Valid types are: 'Fixed', "
335 "'Rotate', or 'Translate'.")
336 }
337
338 m_zones[indx] = zone;
339 zonesElement = zonesElement->NextSiblingElement();
340 }
341}
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:130
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
std::vector< std::string > velocityStr
const NekPoint< NekDouble > origin
std::vector< std::string > displacementStr
std::shared_ptr< ZoneBase > ZoneBaseShPtr
Definition: Zones.h:170
std::shared_ptr< ZoneTranslate > ZoneTranslateShPtr
Definition: Zones.h:313
std::shared_ptr< ZoneRotate > ZoneRotateShPtr
Definition: Zones.h:312
std::shared_ptr< ZoneFixed > ZoneFixedShPtr
Definition: Zones.h:314
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MovementTests::axis, Nektar::MovementTests::displacementStr, Nektar::ParseUtils::GenerateVector(), Nektar::SpatialDomains::MeshGraph::GetDomain(), Nektar::SpatialDomains::MeshGraph::GetSpaceDimension(), m_moveFlag, m_translate, m_zones, Nektar::NekVector< DataType >::Normalize(), Nektar::MovementTests::origin, Nektar::SpatialDomains::ReadTag(), Nektar::MovementTests::velocityStr, and WARNINGL0.

Referenced by Movement().

◆ SetImplicitALEFlag()

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

Definition at line 115 of file Movement.h.

116 {
117 m_ImplicitALESolver = ImplicitALE;
118 }

References m_ImplicitALESolver.

◆ 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 429 of file Movement/Movement.cpp.

430{
431 if (m_zones.size() == 0 && m_interfaces.size() == 0)
432 {
433 return;
434 }
435 TiXmlElement *movement = new TiXmlElement("MOVEMENT");
436 root->LinkEndChild(movement);
437
438 TiXmlElement *zones = new TiXmlElement("ZONES");
439 for (auto &i : m_zones)
440 {
441 const ZoneBaseShPtr z = i.second;
442 const MovementType mtype = z->GetMovementType();
443 std::string label = MovementTypeStr[static_cast<int>(mtype)];
444 TiXmlElement *e = new TiXmlElement(label.substr(0, 1));
445 e->SetAttribute("ID", i.first);
446 std::stringstream s;
447 s << "D[" << z->GetDomainID() << "]";
448 e->SetAttribute("DOMAIN", s.str());
449
450 switch (mtype)
451 {
453 {
454 auto rotate = std::static_pointer_cast<ZoneRotate>(z);
455 e->SetAttribute(
456 "ORIGIN", StripParentheses(rotate->GetOrigin().AsString()));
457 e->SetAttribute("AXIS",
458 StripParentheses(rotate->GetAxis().AsString()));
459 e->SetAttribute("ANGVEL",
460 rotate->GetAngularVelEqn()->GetExpression());
461 }
462 break;
464 {
465 auto translate = std::static_pointer_cast<ZoneTranslate>(z);
466 e->SetAttribute(
467 "XVELOCITY",
468 translate->GetVelocityEquation()[0]->GetExpression());
469 e->SetAttribute(
470 "XDISPLACEMENT",
471 translate->GetDisplacementEquation()[0]->GetExpression());
472 e->SetAttribute(
473 "YVELOCITY",
474 translate->GetVelocityEquation()[1]->GetExpression());
475 e->SetAttribute(
476 "YDISPLACEMENT",
477 translate->GetDisplacementEquation()[1]->GetExpression());
478 e->SetAttribute(
479 "ZVELOCITY",
480 translate->GetVelocityEquation()[2]->GetExpression());
481 e->SetAttribute(
482 "ZDISPLACEMENT",
483 translate->GetDisplacementEquation()[2]->GetExpression());
484 }
485 break;
486 default:
487 break;
488 }
489 zones->LinkEndChild(e);
490 }
491 movement->LinkEndChild(zones);
492
493 TiXmlElement *interfaces = new TiXmlElement("INTERFACES");
494 for (auto &i : m_interfaces)
495 {
496 const std::string interfaceName = i.first.second;
497 TiXmlElement *e = new TiXmlElement("INTERFACE");
498 e->SetAttribute("NAME", interfaceName);
499 const InterfaceShPtr left = i.second->GetLeftInterface();
500 const InterfaceShPtr right = i.second->GetRightInterface();
501 if (left)
502 {
503 TiXmlElement *left_e = new TiXmlElement("L");
504 left_e->SetAttribute("ID", left->GetId());
505 left_e->SetAttribute(
506 "BOUNDARY",
507 "C[" + ParseUtils::GenerateSeqString(left->GetCompositeIDs()) +
508 "]");
509 e->LinkEndChild(left_e);
510 }
511 if (right)
512 {
513 TiXmlElement *right_e = new TiXmlElement("R");
514 right_e->SetAttribute("ID", right->GetId());
515 right_e->SetAttribute(
516 "BOUNDARY",
517 "C[" + ParseUtils::GenerateSeqString(right->GetCompositeIDs()) +
518 "]");
519 e->LinkEndChild(right_e);
520 }
521 interfaces->LinkEndChild(e);
522 }
523 movement->LinkEndChild(interfaces);
524}
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(), CG_Iterations::label, m_interfaces, m_zones, Nektar::SpatialDomains::MovementTypeStr, Nektar::SpatialDomains::StripParentheses(), and Nektar::UnitTests::z().

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 136 of file Movement.h.

Referenced by GetCoordExchangeFlag().

◆ m_DomainBox

Array<OneD, NekDouble> Nektar::SpatialDomains::Movement::m_DomainBox
protected

Definition at line 138 of file Movement.h.

Referenced by DomainBox(), GetDomainBox(), and Movement().

◆ m_DomainLength

Array<OneD, NekDouble> Nektar::SpatialDomains::Movement::m_DomainLength
protected

Definition at line 139 of file Movement.h.

Referenced by DomainBox(), GetDomainLength(), and Movement().

◆ m_ImplicitALESolver

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

Definition at line 135 of file Movement.h.

Referenced by GetImplicitALESolverFlag(), and SetImplicitALEFlag().

◆ m_interfaces

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

◆ m_moved

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

Definition at line 134 of file Movement.h.

Referenced by GetMovedFlag(), and PerformMovement().

◆ m_moveFlag

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

Definition at line 132 of file Movement.h.

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

◆ m_translate

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

Definition at line 133 of file Movement.h.

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

◆ m_zones

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