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

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

Detailed Description

Definition at line 51 of file Movement.h.

Constructor & Destructor Documentation

◆ Movement() [1/2]

Nektar::SpatialDomains::Movement::Movement ( )
inline

Default constructor.

Definition at line 55 of file Movement.h.

56 {
57 }

◆ 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
102 // Don't support interior penalty yet
103 if (pSession->DefinesSolverInfo("DiffusionType"))
104 {
105 ASSERTL0(pSession->GetSolverInfo("DiffusionType") == "LDGNS",
106 "Only LDGNS is supported as the DiffusionType in "
107 "SOLVERINFO when a MOVEMENT block is defined.")
108 }
109 }
110
111 // DEBUG COMMENTS
112 if (movement != nullptr && pSession->DefinesCmdLineArgument("verbose"))
113 {
114 auto comm = pSession->GetComm();
115 if (comm->TreatAsRankZero())
116 {
117 std::cout << "Movement Info:\n";
118 std::cout << "\tNum zones: " << m_zones.size() << "\n";
119 }
120
121 for (auto &zone : m_zones)
122 {
123 auto numEl = zone.second->GetElements().size();
124 comm->GetSpaceComm()->AllReduce(numEl, LibUtilities::ReduceSum);
125
126 // Find shape type if not on this proc
127 int shapeType =
128 !zone.second->GetElements().empty()
129 ? zone.second->GetElements().front()->GetShapeType()
130 : -1;
131 comm->GetSpaceComm()->AllReduce(shapeType, LibUtilities::ReduceMax);
132
133 if (comm->TreatAsRankZero())
134 {
135 std::cout << "\t- " << zone.first << " "
136 << MovementTypeStr[static_cast<int>(
137 zone.second->GetMovementType())]
138 << ": " << numEl << " "
139 << LibUtilities::ShapeTypeMap[shapeType] << "s\n";
140 }
141 }
142
143 if (comm->TreatAsRankZero())
144 {
145 std::cout << "\tNum interfaces: " << m_interfaces.size() << "\n";
146 }
147
148 for (auto &interface : m_interfaces)
149 {
150 auto numLeft =
151 interface.second->GetLeftInterface()->GetEdge().size();
152 auto numRight =
153 interface.second->GetRightInterface()->GetEdge().size();
154 comm->GetSpaceComm()->AllReduce(numLeft, LibUtilities::ReduceSum);
155 comm->GetSpaceComm()->AllReduce(numRight, LibUtilities::ReduceSum);
156
157 // Find shape type if not on this proc
158 int shapeTypeLeft =
159 !interface.second->GetLeftInterface()->GetEdge().empty()
160 ? interface.second->GetLeftInterface()
161 ->GetEdge()
162 .begin()
163 ->second->GetShapeType()
164 : -1;
165 comm->GetSpaceComm()->AllReduce(shapeTypeLeft,
167 int shapeTypeRight =
168 !interface.second->GetRightInterface()->GetEdge().empty()
169 ? interface.second->GetRightInterface()
170 ->GetEdge()
171 .begin()
172 ->second->GetShapeType()
173 : -1;
174 comm->GetSpaceComm()->AllReduce(shapeTypeRight,
176
177 if (comm->TreatAsRankZero())
178 {
179 std::cout << "\t- \"" << interface.first.second << "\": "
180 << interface.second->GetLeftInterface()->GetId()
181 << " (" << numLeft << " "
182 << LibUtilities::ShapeTypeMap[shapeTypeLeft]
183 << "s) <-> "
184 << interface.second->GetRightInterface()->GetId()
185 << " (" << numRight << " "
186 << LibUtilities::ShapeTypeMap[shapeTypeRight]
187 << "s)\n";
188 }
189 }
190
191 comm->GetSpaceComm()->Block();
192 if (comm->TreatAsRankZero())
193 {
194 std::cout << std::endl;
195 }
196 }
197}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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:91
std::map< int, ZoneBaseShPtr > m_zones
Definition: Movement.h:92
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 ASSERTL0, m_interfaces, m_zones, Nektar::SpatialDomains::MovementTypeStr, ReadInterfaces(), ReadZones(), Nektar::LibUtilities::ReduceMax, 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 543 of file Movement/Movement.cpp.

545{
546 m_interfaces[std::make_pair(m_interfaces.size(), name)] =
549}
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 532 of file Movement/Movement.cpp.

533{
534 m_zones[zone->GetId()] = zone;
535 MovementType mtype = zone->GetMovementType();
536 if (mtype != MovementType::eFixed && mtype != MovementType::eNone)
537 {
538 m_moveFlag = true;
539 }
540}
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().

◆ 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().

◆ 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().

◆ PerformMovement()

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

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

506{
507 std::set<int> movedZoneIds;
508 for (auto &zone : m_zones)
509 {
510 if (zone.second->Move(timeStep))
511 {
512 movedZoneIds.insert(zone.first);
513 }
514 }
515
516 // If zone has moved, set all interfaces on that zone to moved.
517 for (auto &interPair : m_interfaces)
518 {
519 int leftId = interPair.second->GetLeftInterface()->GetId();
520 int rightId = interPair.second->GetRightInterface()->GetId();
521
522 if (movedZoneIds.find(leftId) != movedZoneIds.end() ||
523 movedZoneIds.find(rightId) != movedZoneIds.end())
524 {
525 m_zones[leftId]->GetMoved() = true;
526 m_zones[rightId]->GetMoved() = true;
527 }
528 }
529}

References m_interfaces, and m_zones.

◆ ReadInterfaces()

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

Read interfaces given TiXmlDocument.

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

322{
323 ASSERTL0(interfacesTag, "Unable to find INTERFACES tag in file.");
324 TiXmlElement *interfaceElement = interfacesTag->FirstChildElement();
325
326 while (interfaceElement)
327 {
328 ASSERTL0(
329 "INTERFACE" == (std::string)interfaceElement->Value(),
330 "Only INTERFACE tags may be present inside the INTERFACES block.")
331
332 int err;
333
334 std::string name;
335 err = interfaceElement->QueryStringAttribute("NAME", &name);
336 ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface name.");
337 TiXmlElement *sideElement = interfaceElement->FirstChildElement();
338
339 // @TODO: For different interface types have a string attribute type in
340 // @TODO: the INTERFACE element like for NAME above
341
342 Array<OneD, InterfaceShPtr> interfaces(2);
343 std::vector<int> cnt;
344 while (sideElement)
345 {
346 ASSERTL0(
347 cnt.size() < 2,
348 "In INTERFACE NAME " + name +
349 ", only two sides may be present in each INTERFACE block.")
350
351 int indx;
352 err = sideElement->QueryIntAttribute("ID", &indx);
353 ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface ID.");
354
355 std::string boundaryStr;
356 int boundaryErr =
357 sideElement->QueryStringAttribute("BOUNDARY", &boundaryStr);
358
359 CompositeMap boundaryEdge;
360 std::string indxStr;
361 if (boundaryErr == TIXML_SUCCESS)
362 {
363 indxStr = ReadTag(boundaryStr);
364 meshGraph->GetCompositeList(indxStr, boundaryEdge);
365 }
366
367 // Sets location in interface pair to 0 for left, and 1 for right
368 auto sideElVal = sideElement->ValueStr();
369 if (sideElVal == "LEFT" || sideElVal == "L")
370 {
371 cnt.emplace_back(0);
372 }
373 else if (sideElVal == "RIGHT" || sideElVal == "R")
374 {
375 cnt.emplace_back(1);
376 }
377 else
378 {
380 sideElement->ValueStr() +
381 " is not a valid interface side for interface "
382 "NAME " +
383 name + ". Please only use LEFT or RIGHT.")
384 }
385
386 interfaces[cnt[cnt.size() - 1]] =
388 indx, boundaryEdge));
389
390 sideElement = sideElement->NextSiblingElement();
391 }
392
393 ASSERTL0(std::accumulate(cnt.begin(), cnt.end(), 0) == 1,
394 "You must have only one LEFT and one RIGHT side"
395 " present in interface NAME " +
396 name)
397
398 m_interfaces[std::make_pair(m_interfaces.size(), name)] =
400 interfaces[0], interfaces[1]));
401 interfaceElement = interfaceElement->NextSiblingElement();
402 }
403}
#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

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

201{
202 int coordDim = meshGraph->GetSpaceDimension();
203
204 ASSERTL0(zonesTag, "Unable to find ZONES tag in file.");
205 TiXmlElement *zonesElement = zonesTag->FirstChildElement();
206 while (zonesElement)
207 {
208 std::string zoneType = zonesElement->Value();
209
210 int err;
211 int indx;
212
213 err = zonesElement->QueryIntAttribute("ID", &indx);
214 ASSERTL0(err == TIXML_SUCCESS, "Unable to read zone ID.");
215
216 std::string interfaceDomainStr;
217 err = zonesElement->QueryStringAttribute("DOMAIN", &interfaceDomainStr);
218 ASSERTL0(err == TIXML_SUCCESS, "Unable to read zone domain.");
219
220 auto &domains = meshGraph->GetDomain();
221 auto domFind = std::stoi(ReadTag(interfaceDomainStr));
222 std::map<int, CompositeSharedPtr> domain;
223 if (domains.find(domFind) != domains.end())
224 {
225 domain = domains.at(domFind);
226 }
227
228 ZoneBaseShPtr zone;
229 if (zoneType == "F" || zoneType == "FIXED")
230 {
232 indx, domFind, domain, coordDim));
233 }
234 else if (zoneType == "R" || zoneType == "ROTATE" ||
235 zoneType == "ROTATING")
236 {
237 std::string originStr;
238 err = zonesElement->QueryStringAttribute("ORIGIN", &originStr);
239 ASSERTL0(err == TIXML_SUCCESS, "Unable to read origin.");
240 std::vector<NekDouble> originVec;
241 ParseUtils::GenerateVector(originStr, originVec);
242 NekPoint<NekDouble> origin =
243 NekPoint<NekDouble>(originVec[0], originVec[1], originVec[2]);
244
245 std::string axisStr;
246 err = zonesElement->QueryStringAttribute("AXIS", &axisStr);
247 ASSERTL0(err == TIXML_SUCCESS, "Unable to read axis.");
248 DNekVec axis(axisStr);
249 axis.Normalize();
250
251 std::string angularVelStr;
252 err = zonesElement->QueryStringAttribute("ANGVEL", &angularVelStr);
253 ASSERTL0(err == TIXML_SUCCESS, "Unable to read angular velocity.");
254
255 LibUtilities::EquationSharedPtr angularVelEqn =
257 pSession->GetInterpreter(), angularVelStr);
258
260 indx, domFind, domain, coordDim, origin, axis, angularVelEqn));
261
262 m_moveFlag = true;
263 }
264 else if (zoneType == "T" || zoneType == "TRANSLATE" ||
265 zoneType == "TRANSLATING")
266 {
267 std::string velocityStr;
268 err = zonesElement->QueryStringAttribute("VELOCITY", &velocityStr);
269 ASSERTL0(err == TIXML_SUCCESS, "Unable to read direction.");
270 std::vector<NekDouble> velocity;
272
273 zone = ZoneTranslateShPtr(
275 indx, domFind, domain, coordDim, velocity));
276
277 m_moveFlag = true;
278 }
279 else if (zoneType == "P" || zoneType == "PRESCRIBED")
280 {
281 std::string xDeformStr;
282 err = zonesElement->QueryStringAttribute("XDEFORM", &xDeformStr);
283 ASSERTL0(err == TIXML_SUCCESS, "Unable to read x deform equation.");
285 std::make_shared<LibUtilities::Equation>(
286 pSession->GetInterpreter(), xDeformStr);
287
288 std::string yDeformStr;
289 err = zonesElement->QueryStringAttribute("YDEFORM", &yDeformStr);
290 ASSERTL0(err == TIXML_SUCCESS, "Unable to read y deform equation.");
292 std::make_shared<LibUtilities::Equation>(
293 pSession->GetInterpreter(), yDeformStr);
294
295 std::string zDeformStr;
296 err = zonesElement->QueryStringAttribute("ZDEFORM", &zDeformStr);
297 ASSERTL0(err == TIXML_SUCCESS, "Unable to read z deform equation.");
299 std::make_shared<LibUtilities::Equation>(
300 pSession->GetInterpreter(), zDeformStr);
301
302 zone = ZonePrescribeShPtr(
304 indx, domFind, domain, coordDim, xDeformEqn, yDeformEqn,
305 zDeformEqn));
306
307 m_moveFlag = true;
308 }
309 else
310 {
311 WARNINGL0(false, "Zone type '" + zoneType +
312 "' is unsupported. Valid types are: 'Fixed', "
313 "'Rotate', 'Translate', or 'Prescribe'.")
314 }
315
316 m_zones[indx] = zone;
317 zonesElement = zonesElement->NextSiblingElement();
318 }
319}
#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
const std::vector< NekDouble > velocity
const NekPoint< NekDouble > origin
std::shared_ptr< ZoneBase > ZoneBaseShPtr
Definition: Zones.h:144
std::shared_ptr< ZonePrescribe > ZonePrescribeShPtr
Definition: Zones.h:364
std::shared_ptr< ZoneTranslate > ZoneTranslateShPtr
Definition: Zones.h:363
std::shared_ptr< ZoneRotate > ZoneRotateShPtr
Definition: Zones.h:362
std::shared_ptr< ZoneFixed > ZoneFixedShPtr
Definition: Zones.h:365
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48

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

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

407{
408 if (m_zones.size() == 0 && m_interfaces.size() == 0)
409 {
410 return;
411 }
412 TiXmlElement *movement = new TiXmlElement("MOVEMENT");
413 root->LinkEndChild(movement);
414
415 TiXmlElement *zones = new TiXmlElement("ZONES");
416 for (auto &i : m_zones)
417 {
418 const ZoneBaseShPtr z = i.second;
419 const MovementType mtype = z->GetMovementType();
420 std::string label = MovementTypeStr[static_cast<int>(mtype)];
421 TiXmlElement *e = new TiXmlElement(label.substr(0, 1));
422 e->SetAttribute("ID", i.first);
423 std::stringstream s;
424 s << "D[" << z->GetDomainID() << "]";
425 e->SetAttribute("DOMAIN", s.str());
426
427 switch (mtype)
428 {
430 {
431 auto rotate = std::static_pointer_cast<ZoneRotate>(z);
432 e->SetAttribute(
433 "ORIGIN", StripParentheses(rotate->GetOrigin().AsString()));
434 e->SetAttribute("AXIS",
435 StripParentheses(rotate->GetAxis().AsString()));
436 e->SetAttribute("ANGVEL",
437 rotate->GetAngularVelEqn()->GetExpression());
438 }
439 break;
441 {
442 auto translate = std::static_pointer_cast<ZoneTranslate>(z);
443 const std::vector<NekDouble> vel = translate->GetVel();
444 std::stringstream vel_s;
445 vel_s << vel[0] << ", " << vel[1] << ", " << vel[2];
446 e->SetAttribute("VELOCITY", vel_s.str());
447 }
448 break;
450 {
451 auto prescribe = std::static_pointer_cast<ZonePrescribe>(z);
452 e->SetAttribute(
453 "XDEFORM",
454 prescribe->GetXDeformEquation()->GetExpression());
455 e->SetAttribute(
456 "YDEFORM",
457 prescribe->GetYDeformEquation()->GetExpression());
458 e->SetAttribute(
459 "ZDEFORM",
460 prescribe->GetZDeformEquation()->GetExpression());
461 }
462 break;
463 default:
464 break;
465 }
466 zones->LinkEndChild(e);
467 }
468 movement->LinkEndChild(zones);
469
470 TiXmlElement *interfaces = new TiXmlElement("INTERFACES");
471 for (auto &i : m_interfaces)
472 {
473 const std::string interfaceName = i.first.second;
474 TiXmlElement *e = new TiXmlElement("INTERFACE");
475 e->SetAttribute("NAME", interfaceName);
476 const InterfaceShPtr left = i.second->GetLeftInterface();
477 const InterfaceShPtr right = i.second->GetRightInterface();
478 if (left)
479 {
480 TiXmlElement *left_e = new TiXmlElement("L");
481 left_e->SetAttribute("ID", left->GetId());
482 left_e->SetAttribute(
483 "BOUNDARY",
484 "C[" + ParseUtils::GenerateSeqString(left->GetCompositeIDs()) +
485 "]");
486 e->LinkEndChild(left_e);
487 }
488 if (right)
489 {
490 TiXmlElement *right_e = new TiXmlElement("R");
491 right_e->SetAttribute("ID", right->GetId());
492 right_e->SetAttribute(
493 "BOUNDARY",
494 "C[" + ParseUtils::GenerateSeqString(right->GetCompositeIDs()) +
495 "]");
496 e->LinkEndChild(right_e);
497 }
498 interfaces->LinkEndChild(e);
499 }
500 movement->LinkEndChild(interfaces);
501}
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::ePrescribe, 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_interfaces

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

◆ m_moveFlag

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

Definition at line 93 of file Movement.h.

Referenced by AddZone(), and ReadZones().

◆ m_zones

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

Definition at line 92 of file Movement.h.

Referenced by AddZone(), GetZones(), Movement(), PerformMovement(), ReadZones(), and WriteMovement().