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 (const LibUtilities::SessionReaderSharedPtr &pSession, MeshGraph *meshGraph)
 Constructor. More...
 
 ~Movement ()=default
 Default destructor. More...
 
const InterfaceCollectionGetInterfaces () const
 
const std::map< int, ZoneBaseShPtr > & GetZones () const
 
void PerformMovement (NekDouble timeStep)
 

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

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

Constructor.

Definition at line 62 of file Movement.cpp.

64 {
65  TiXmlNode *nektar = pSession->GetElement("NEKTAR");
66  if (nektar == nullptr)
67  {
68  return;
69  }
70 
71  TiXmlNode *movement = nektar->FirstChild("MOVEMENT");
72  if (movement != nullptr)
73  {
74  bool zones = movement->FirstChild("ZONES") != nullptr;
75  if (zones)
76  {
77  ReadZones(pSession->GetElement("NEKTAR/MOVEMENT/ZONES"), meshGraph,
78  pSession);
79  }
80 
81  bool interfaces = movement->FirstChild("INTERFACES") != nullptr;
82  if (interfaces)
83  {
84  ReadInterfaces(pSession->GetElement("NEKTAR/MOVEMENT/INTERFACES"),
85  meshGraph);
86  }
87 
88  ASSERTL0(zones == interfaces,
89  "Only one of ZONES or INTERFACES present in the MOVEMENT "
90  "block.")
91 
92  // Don't support interior penalty yet
93  if (pSession->DefinesSolverInfo("DiffusionType"))
94  {
95  ASSERTL0(pSession->GetSolverInfo("DiffusionType") == "LDGNS",
96  "Only LDGNS is supported as the DiffusionType in "
97  "SOLVERINFO when a MOVEMENT block is defined.")
98  }
99  }
100 
101  // DEBUG COMMENTS
102  if (movement != nullptr && pSession->DefinesCmdLineArgument("verbose"))
103  {
104  auto comm = pSession->GetComm();
105  if (comm->TreatAsRankZero())
106  {
107  std::cout << "Movement Info:\n";
108  std::cout << "\tNum zones: " << m_zones.size() << "\n";
109  }
110 
111  for (auto &zone : m_zones)
112  {
113  auto numEl = zone.second->GetElements().size();
114  comm->AllReduce(numEl, LibUtilities::ReduceSum);
115 
116  // Find shape type if not on this proc
117  int shapeType =
118  !zone.second->GetElements().empty()
119  ? zone.second->GetElements().front()->GetShapeType()
120  : -1;
121  comm->AllReduce(shapeType, LibUtilities::ReduceMax);
122 
123  if (comm->TreatAsRankZero())
124  {
125  std::cout << "\t- " << zone.first << " "
126  << MovementTypeStr[static_cast<int>(
127  zone.second->GetMovementType())]
128  << ": " << numEl << " "
129  << LibUtilities::ShapeTypeMap[shapeType] << "s\n";
130  }
131  }
132 
133  if (comm->TreatAsRankZero())
134  {
135  std::cout << "\tNum interfaces: " << m_interfaces.size() << "\n";
136  }
137 
138  for (auto &interface : m_interfaces)
139  {
140  auto numLeft =
141  interface.second->GetLeftInterface()->GetEdge().size();
142  auto numRight =
143  interface.second->GetRightInterface()->GetEdge().size();
144  comm->AllReduce(numLeft, LibUtilities::ReduceSum);
145  comm->AllReduce(numRight, LibUtilities::ReduceSum);
146 
147  // Find shape type if not on this proc
148  int shapeTypeLeft =
149  !interface.second->GetLeftInterface()->GetEdge().empty()
150  ? interface.second->GetLeftInterface()
151  ->GetEdge()
152  .begin()
153  ->second->GetShapeType()
154  : -1;
155  comm->AllReduce(shapeTypeLeft, LibUtilities::ReduceMax);
156  int shapeTypeRight =
157  !interface.second->GetRightInterface()->GetEdge().empty()
158  ? interface.second->GetRightInterface()
159  ->GetEdge()
160  .begin()
161  ->second->GetShapeType()
162  : -1;
163  comm->AllReduce(shapeTypeRight, LibUtilities::ReduceMax);
164 
165  if (comm->TreatAsRankZero())
166  {
167  std::cout << "\t- \"" << interface.first.second << "\": "
168  << interface.second->GetLeftInterface()->GetId()
169  << " (" << numLeft << " "
170  << LibUtilities::ShapeTypeMap[shapeTypeLeft]
171  << "s) <-> "
172  << interface.second->GetRightInterface()->GetId()
173  << " (" << numRight << " "
174  << LibUtilities::ShapeTypeMap[shapeTypeRight]
175  << "s)\n";
176  }
177  }
178 
179  comm->Block();
180  if (comm->TreatAsRankZero())
181  {
182  std::cout << std::endl;
183  }
184  }
185 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
void ReadZones(TiXmlElement *zonesTag, MeshGraph *meshGraph, const LibUtilities::SessionReaderSharedPtr &pSession)
Read zones given TiXmlDocument.
Definition: Movement.cpp:187
void ReadInterfaces(TiXmlElement *interfacesTag, MeshGraph *meshGraph)
Read interfaces given TiXmlDocument.
Definition: Movement.cpp:310
InterfaceCollection m_interfaces
Definition: Movement.h:76
std::map< int, ZoneBaseShPtr > m_zones
Definition: Movement.h:77
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:79
const std::string MovementTypeStr[]
Map of zone movement type to movement type string.
Definition: Zones.h:58

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

◆ GetInterfaces()

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

Definition at line 62 of file Movement.h.

63  {
64  return m_interfaces;
65  }

References m_interfaces.

◆ GetZones()

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

Definition at line 67 of file Movement.h.

68  {
69  return m_zones;
70  }

References m_zones.

◆ PerformMovement()

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

Definition at line 395 of file Movement.cpp.

396 {
397  std::set<int> movedZoneIds;
398  for (auto &zone : m_zones)
399  {
400  if (zone.second->Move(timeStep))
401  {
402  movedZoneIds.insert(zone.first);
403  }
404  }
405 
406  // If zone has moved, set all interfaces on that zone to moved.
407  for (auto &interPair : m_interfaces)
408  {
409  int leftId = interPair.second->GetLeftInterface()->GetId();
410  int rightId = interPair.second->GetRightInterface()->GetId();
411 
412  if (movedZoneIds.find(leftId) != movedZoneIds.end() ||
413  movedZoneIds.find(rightId) != movedZoneIds.end())
414  {
415  m_zones[leftId]->GetMoved() = true;
416  m_zones[rightId]->GetMoved() = true;
417  }
418  }
419 }

References m_interfaces, and m_zones.

◆ ReadInterfaces()

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

Read interfaces given TiXmlDocument.

Definition at line 310 of file Movement.cpp.

311 {
312  ASSERTL0(interfacesTag, "Unable to find INTERFACES tag in file.");
313  TiXmlElement *interfaceElement = interfacesTag->FirstChildElement();
314 
315  while (interfaceElement)
316  {
317  ASSERTL0(
318  "INTERFACE" == (std::string)interfaceElement->Value(),
319  "Only INTERFACE tags may be present inside the INTERFACES block.")
320 
321  int err;
322 
323  std::string name;
324  err = interfaceElement->QueryStringAttribute("NAME", &name);
325  ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface name.");
326  TiXmlElement *sideElement = interfaceElement->FirstChildElement();
327 
328  // @TODO: For different interface types have a string attribute type in
329  // @TODO: the INTERFACE element like for NAME above
330 
331  Array<OneD, InterfaceShPtr> interfaces(2);
332  std::vector<int> cnt;
333  while (sideElement)
334  {
335  ASSERTL0(
336  cnt.size() < 2,
337  "In INTERFACE NAME " + name +
338  ", only two sides may be present in each INTERFACE block.")
339 
340  int indx;
341  err = sideElement->QueryIntAttribute("ID", &indx);
342  ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface ID.");
343 
344  std::string boundaryStr;
345  int boundaryErr =
346  sideElement->QueryStringAttribute("BOUNDARY", &boundaryStr);
347 
348  CompositeMap boundaryEdge;
349  if (boundaryErr == TIXML_SUCCESS)
350  {
351  std::string indxStr = ReadTag(boundaryStr);
352  meshGraph->GetCompositeList(indxStr, boundaryEdge);
353  }
354 
355  // Sets location in interface pair to 0 for left, and 1 for right
356  auto sideElVal = sideElement->ValueStr();
357  if (sideElVal == "LEFT" || sideElVal == "L")
358  {
359  cnt.emplace_back(0);
360  }
361  else if (sideElVal == "RIGHT" || sideElVal == "R")
362  {
363  cnt.emplace_back(1);
364  }
365  else
366  {
368  sideElement->ValueStr() +
369  " is not a valid interface side for interface "
370  "NAME " +
371  name + ". Please only use LEFT or RIGHT.")
372  }
373 
374  interfaces[cnt[cnt.size() - 1]] =
376  indx, boundaryEdge));
377 
378  sideElement = sideElement->NextSiblingElement();
379  }
380 
381  ASSERTL0(std::accumulate(cnt.begin(), cnt.end(), 0) == 1,
382  "You must have only one LEFT and one RIGHT side"
383  " present in interface NAME " +
384  name)
385 
386  m_interfaces[std::make_pair(m_interfaces.size(), name)] =
388  interfaces[0], interfaces[1]));
389  interfaceElement = interfaceElement->NextSiblingElement();
390  }
391 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< InterfacePair > InterfacePairShPtr
std::shared_ptr< Interface > InterfaceShPtr
static std::string ReadTag(std::string &tagStr)
Definition: Movement.cpp:47
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:138

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

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

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

Referenced by Movement().

Member Data Documentation

◆ m_interfaces

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

Definition at line 76 of file Movement.h.

Referenced by GetInterfaces(), Movement(), PerformMovement(), and ReadInterfaces().

◆ m_moveFlag

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

Definition at line 78 of file Movement.h.

Referenced by ReadZones().

◆ m_zones

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

Definition at line 77 of file Movement.h.

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