Nektar++
Movement.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: Movement.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: This file contains the base class for implementing
33// non-conformal geometry using the Movement object
34//
35////////////////////////////////////////////////////////////////////////////////
36
40#include <tinyxml.h>
41
42namespace Nektar
43{
44namespace SpatialDomains
45{
46
47std::string static inline ReadTag(std::string &tagStr)
48{
49 std::string::size_type indxBeg = tagStr.find_first_of('[') + 1;
50 std::string::size_type indxEnd = tagStr.find_last_of(']') - 1;
51
53 indxBeg <= indxEnd,
54 (std::string("Error reading interface region definition:") + tagStr)
55 .c_str());
56
57 std::string indxStr = tagStr.substr(indxBeg, indxEnd - indxBeg + 1);
58
59 return indxStr;
60}
61
63 MeshGraph *meshGraph)
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 TiXmlElement *zonesTag =
78 pSession->GetElement("NEKTAR/MOVEMENT/ZONES");
79 ReadZones(zonesTag, meshGraph, pSession);
80 }
81
82 bool interfaces = movement->FirstChild("INTERFACES") != nullptr;
83 if (interfaces)
84 {
85 TiXmlElement *interfacesTag =
86 pSession->GetElement("NEKTAR/MOVEMENT/INTERFACES");
87 ReadInterfaces(interfacesTag, meshGraph);
88 }
89
90 ASSERTL0(zones == interfaces,
91 "Only one of ZONES or INTERFACES present in the MOVEMENT "
92 "block.")
93
94 // Don't support interior penalty yet
95 if (pSession->DefinesSolverInfo("DiffusionType"))
96 {
97 ASSERTL0(pSession->GetSolverInfo("DiffusionType") == "LDGNS",
98 "Only LDGNS is supported as the DiffusionType in "
99 "SOLVERINFO when a MOVEMENT block is defined.")
100 }
101 }
102
103 // DEBUG COMMENTS
104 if (movement != nullptr && pSession->DefinesCmdLineArgument("verbose"))
105 {
106 auto comm = pSession->GetComm();
107 if (comm->TreatAsRankZero())
108 {
109 std::cout << "Movement Info:\n";
110 std::cout << "\tNum zones: " << m_zones.size() << "\n";
111 }
112
113 for (auto &zone : m_zones)
114 {
115 auto numEl = zone.second->GetElements().size();
116 comm->GetSpaceComm()->AllReduce(numEl, LibUtilities::ReduceSum);
117
118 // Find shape type if not on this proc
119 int shapeType =
120 !zone.second->GetElements().empty()
121 ? zone.second->GetElements().front()->GetShapeType()
122 : -1;
123 comm->GetSpaceComm()->AllReduce(shapeType, LibUtilities::ReduceMax);
124
125 if (comm->TreatAsRankZero())
126 {
127 std::cout << "\t- " << zone.first << " "
128 << MovementTypeStr[static_cast<int>(
129 zone.second->GetMovementType())]
130 << ": " << numEl << " "
131 << LibUtilities::ShapeTypeMap[shapeType] << "s\n";
132 }
133 }
134
135 if (comm->TreatAsRankZero())
136 {
137 std::cout << "\tNum interfaces: " << m_interfaces.size() << "\n";
138 }
139
140 for (auto &interface : m_interfaces)
141 {
142 auto numLeft =
143 interface.second->GetLeftInterface()->GetEdge().size();
144 auto numRight =
145 interface.second->GetRightInterface()->GetEdge().size();
146 comm->GetSpaceComm()->AllReduce(numLeft, LibUtilities::ReduceSum);
147 comm->GetSpaceComm()->AllReduce(numRight, LibUtilities::ReduceSum);
148
149 // Find shape type if not on this proc
150 int shapeTypeLeft =
151 !interface.second->GetLeftInterface()->GetEdge().empty()
152 ? interface.second->GetLeftInterface()
153 ->GetEdge()
154 .begin()
155 ->second->GetShapeType()
156 : -1;
157 comm->GetSpaceComm()->AllReduce(shapeTypeLeft,
159 int shapeTypeRight =
160 !interface.second->GetRightInterface()->GetEdge().empty()
161 ? interface.second->GetRightInterface()
162 ->GetEdge()
163 .begin()
164 ->second->GetShapeType()
165 : -1;
166 comm->GetSpaceComm()->AllReduce(shapeTypeRight,
168
169 if (comm->TreatAsRankZero())
170 {
171 std::cout << "\t- \"" << interface.first.second << "\": "
172 << interface.second->GetLeftInterface()->GetId()
173 << " (" << numLeft << " "
174 << LibUtilities::ShapeTypeMap[shapeTypeLeft]
175 << "s) <-> "
176 << interface.second->GetRightInterface()->GetId()
177 << " (" << numRight << " "
178 << LibUtilities::ShapeTypeMap[shapeTypeRight]
179 << "s)\n";
180 }
181 }
182
183 comm->GetSpaceComm()->Block();
184 if (comm->TreatAsRankZero())
185 {
186 std::cout << std::endl;
187 }
188 }
189}
190
191void Movement::ReadZones(TiXmlElement *zonesTag, MeshGraph *meshGraph,
193{
194 int coordDim = meshGraph->GetSpaceDimension();
195
196 ASSERTL0(zonesTag, "Unable to find ZONES tag in file.");
197 TiXmlElement *zonesElement = zonesTag->FirstChildElement();
198 while (zonesElement)
199 {
200 std::string zoneType = zonesElement->Value();
201
202 int err;
203 int indx;
204
205 err = zonesElement->QueryIntAttribute("ID", &indx);
206 ASSERTL0(err == TIXML_SUCCESS, "Unable to read zone ID.");
207
208 std::string interfaceDomainStr;
209 err = zonesElement->QueryStringAttribute("DOMAIN", &interfaceDomainStr);
210 ASSERTL0(err == TIXML_SUCCESS, "Unable to read zone domain.");
211
212 auto &domains = meshGraph->GetDomain();
213 auto domFind = stoi(ReadTag(interfaceDomainStr));
214 std::map<int, CompositeSharedPtr> domain;
215 if (domains.find(domFind) != domains.end())
216 {
217 domain = domains.at(domFind);
218 }
219
220 ZoneBaseShPtr zone;
221 if (zoneType == "F" || zoneType == "FIXED")
222 {
224 indx, domain, coordDim));
225 }
226 else if (zoneType == "R" || zoneType == "ROTATE" ||
227 zoneType == "ROTATING")
228 {
229 std::string originStr;
230 err = zonesElement->QueryStringAttribute("ORIGIN", &originStr);
231 ASSERTL0(err == TIXML_SUCCESS, "Unable to read origin.");
232 std::vector<NekDouble> originVec;
233 ParseUtils::GenerateVector(originStr, originVec);
234 NekPoint<NekDouble> origin =
235 NekPoint<NekDouble>(originVec[0], originVec[1], originVec[2]);
236
237 std::string axisStr;
238 err = zonesElement->QueryStringAttribute("AXIS", &axisStr);
239 ASSERTL0(err == TIXML_SUCCESS, "Unable to read axis.");
240 DNekVec axis(axisStr);
241 axis.Normalize();
242
243 std::string angularVelStr;
244 err = zonesElement->QueryStringAttribute("ANGVEL", &angularVelStr);
245 ASSERTL0(err == TIXML_SUCCESS, "Unable to read angular velocity.");
246
247 LibUtilities::EquationSharedPtr angularVelEqn =
249 pSession->GetInterpreter(), angularVelStr);
250
252 indx, domain, coordDim, origin, axis, angularVelEqn));
253
254 m_moveFlag = true;
255 }
256 else if (zoneType == "T" || zoneType == "TRANSLATE" ||
257 zoneType == "TRANSLATING")
258 {
259 std::string velocityStr;
260 err = zonesElement->QueryStringAttribute("VELOCITY", &velocityStr);
261 ASSERTL0(err == TIXML_SUCCESS, "Unable to read direction.");
262 std::vector<NekDouble> velocity;
263 ParseUtils::GenerateVector(velocityStr, velocity);
264
265 zone = ZoneTranslateShPtr(
267 indx, domain, coordDim, velocity));
268
269 m_moveFlag = true;
270 }
271 else if (zoneType == "P" || zoneType == "PRESCRIBED")
272 {
273 std::string xDeformStr;
274 err = zonesElement->QueryStringAttribute("XDEFORM", &xDeformStr);
275 ASSERTL0(err == TIXML_SUCCESS, "Unable to read x deform equation.");
277 std::make_shared<LibUtilities::Equation>(
278 pSession->GetInterpreter(), xDeformStr);
279
280 std::string yDeformStr;
281 err = zonesElement->QueryStringAttribute("YDEFORM", &yDeformStr);
282 ASSERTL0(err == TIXML_SUCCESS, "Unable to read y deform equation.");
284 std::make_shared<LibUtilities::Equation>(
285 pSession->GetInterpreter(), yDeformStr);
286
287 std::string zDeformStr;
288 err = zonesElement->QueryStringAttribute("ZDEFORM", &zDeformStr);
289 ASSERTL0(err == TIXML_SUCCESS, "Unable to read z deform equation.");
291 std::make_shared<LibUtilities::Equation>(
292 pSession->GetInterpreter(), zDeformStr);
293
294 zone = ZonePrescribeShPtr(
296 indx, domain, coordDim, xDeformEqn, yDeformEqn,
297 zDeformEqn));
298
299 m_moveFlag = true;
300 }
301 else
302 {
303 WARNINGL0(false, "Zone type '" + zoneType +
304 "' is unsupported. Valid types are: 'Fixed', "
305 "'Rotate', 'Translate', or 'Prescribe'.")
306 }
307
308 m_zones[indx] = zone;
309 zonesElement = zonesElement->NextSiblingElement();
310 }
311}
312
313void Movement::ReadInterfaces(TiXmlElement *interfacesTag, MeshGraph *meshGraph)
314{
315 ASSERTL0(interfacesTag, "Unable to find INTERFACES tag in file.");
316 TiXmlElement *interfaceElement = interfacesTag->FirstChildElement();
317
318 while (interfaceElement)
319 {
320 ASSERTL0(
321 "INTERFACE" == (std::string)interfaceElement->Value(),
322 "Only INTERFACE tags may be present inside the INTERFACES block.")
323
324 int err;
325
326 std::string name;
327 err = interfaceElement->QueryStringAttribute("NAME", &name);
328 ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface name.");
329 TiXmlElement *sideElement = interfaceElement->FirstChildElement();
330
331 // @TODO: For different interface types have a string attribute type in
332 // @TODO: the INTERFACE element like for NAME above
333
334 Array<OneD, InterfaceShPtr> interfaces(2);
335 std::vector<int> cnt;
336 while (sideElement)
337 {
338 ASSERTL0(
339 cnt.size() < 2,
340 "In INTERFACE NAME " + name +
341 ", only two sides may be present in each INTERFACE block.")
342
343 int indx;
344 err = sideElement->QueryIntAttribute("ID", &indx);
345 ASSERTL0(err == TIXML_SUCCESS, "Unable to read interface ID.");
346
347 std::string boundaryStr;
348 int boundaryErr =
349 sideElement->QueryStringAttribute("BOUNDARY", &boundaryStr);
350
351 CompositeMap boundaryEdge;
352 if (boundaryErr == TIXML_SUCCESS)
353 {
354 std::string indxStr = ReadTag(boundaryStr);
355 meshGraph->GetCompositeList(indxStr, boundaryEdge);
356 }
357
358 // Sets location in interface pair to 0 for left, and 1 for right
359 auto sideElVal = sideElement->ValueStr();
360 if (sideElVal == "LEFT" || sideElVal == "L")
361 {
362 cnt.emplace_back(0);
363 }
364 else if (sideElVal == "RIGHT" || sideElVal == "R")
365 {
366 cnt.emplace_back(1);
367 }
368 else
369 {
371 sideElement->ValueStr() +
372 " is not a valid interface side for interface "
373 "NAME " +
374 name + ". Please only use LEFT or RIGHT.")
375 }
376
377 interfaces[cnt[cnt.size() - 1]] =
379 indx, boundaryEdge));
380
381 sideElement = sideElement->NextSiblingElement();
382 }
383
384 ASSERTL0(std::accumulate(cnt.begin(), cnt.end(), 0) == 1,
385 "You must have only one LEFT and one RIGHT side"
386 " present in interface NAME " +
387 name)
388
389 m_interfaces[std::make_pair(m_interfaces.size(), name)] =
391 interfaces[0], interfaces[1]));
392 interfaceElement = interfaceElement->NextSiblingElement();
393 }
394}
395
396// Acts as a placeholder for when ALE function and moving geometry capability
397// is added. Currently unused.
399{
400 std::set<int> movedZoneIds;
401 for (auto &zone : m_zones)
402 {
403 if (zone.second->Move(timeStep))
404 {
405 movedZoneIds.insert(zone.first);
406 }
407 }
408
409 // If zone has moved, set all interfaces on that zone to moved.
410 for (auto &interPair : m_interfaces)
411 {
412 int leftId = interPair.second->GetLeftInterface()->GetId();
413 int rightId = interPair.second->GetRightInterface()->GetId();
414
415 if (movedZoneIds.find(leftId) != movedZoneIds.end() ||
416 movedZoneIds.find(rightId) != movedZoneIds.end())
417 {
418 m_zones[leftId]->GetMoved() = true;
419 m_zones[rightId]->GetMoved() = true;
420 }
421 }
422}
423
424} // namespace SpatialDomains
425} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:222
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:131
Base class for a spectral/hp element mesh.
Definition: MeshGraph.h:183
std::map< int, std::map< int, CompositeSharedPtr > > & GetDomain()
Definition: MeshGraph.h:271
void GetCompositeList(const std::string &compositeStr, CompositeMap &compositeVector) const
Definition: MeshGraph.cpp:601
int GetSpaceDimension()
Dimension of the space (can be a 1D curve in 3D space).
Definition: MeshGraph.h:227
void PerformMovement(NekDouble timeStep)
Definition: Movement.cpp:398
Movement(const LibUtilities::SessionReaderSharedPtr &pSession, MeshGraph *meshGraph)
Constructor.
Definition: Movement.cpp:62
void ReadZones(TiXmlElement *zonesTag, MeshGraph *meshGraph, const LibUtilities::SessionReaderSharedPtr &pSession)
Read zones given TiXmlDocument.
Definition: Movement.cpp:191
void ReadInterfaces(TiXmlElement *interfacesTag, MeshGraph *meshGraph)
Read interfaces given TiXmlDocument.
Definition: Movement.cpp:313
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
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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:313
std::shared_ptr< ZoneTranslate > ZoneTranslateShPtr
Definition: Zones.h:312
std::shared_ptr< InterfacePair > InterfacePairShPtr
std::shared_ptr< ZoneRotate > ZoneRotateShPtr
Definition: Zones.h:311
std::shared_ptr< Interface > InterfaceShPtr
std::shared_ptr< ZoneFixed > ZoneFixedShPtr
Definition: Zones.h:314
static std::string ReadTag(std::string &tagStr)
Definition: Movement.cpp:47
const std::string MovementTypeStr[]
Map of zone movement type to movement type string.
Definition: Zones.h:58
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:138
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble