Nektar++
Functions | Variables
Nektar::MovementTests Namespace Reference

Functions

LibUtilities::SessionReaderSharedPtr CreateSession ()
 
SpatialDomains::ZoneBaseShPtr CreateZone (SpatialDomains::MovementType type, int zoneID, int domainID, LibUtilities::InterpreterSharedPtr interpreter)
 Produce dummy Zone objects, containing empty domain pointers. More...
 
SpatialDomains::InterfaceShPtr CreateInterface (int interfaceID, std::vector< int > compositeIDs)
 Produce dummy Interface objects, containing empty domain pointers. More...
 
 BOOST_AUTO_TEST_CASE (TestAddGetZones)
 
 BOOST_AUTO_TEST_CASE (TestAddGetInterfaces)
 
 BOOST_AUTO_TEST_CASE (TestWriteMovement)
 

Variables

const std::string angVelStr = "0.1*t"
 
const std::string xEqnStr = "0.1*x - 0.1*t"
 
const std::string yEqnStr = "0.1*y^2 - x"
 
const std::string zEqnStr = "sqrt(t)"
 
std::vector< std::string > velocityStr = {"1.0", "2.0", "3.0"}
 
std::vector< std::string > displacementStr = {"1.0", "2.0", "3.0"}
 
const NekPoint< NekDoubleorigin = {1., 2., 3.}
 
const DNekVec axis = {1., 2., 3.}
 

Function Documentation

◆ BOOST_AUTO_TEST_CASE() [1/3]

Nektar::MovementTests::BOOST_AUTO_TEST_CASE ( TestAddGetInterfaces  )

Definition at line 147 of file TestMovement.cpp.

148{
151 interface2 = CreateInterface(1, {1, 2, 3}),
152 interface3 = CreateInterface(3, {4});
153 m.AddInterface("north", interface1, interface2);
154 m.AddInterface("east", interface2, interface3);
155 m.AddInterface("south", interface3, interface1);
157 BOOST_TEST(interfaces.size() == 3);
158 SpatialDomains::InterfacePairShPtr north = interfaces.at(
159 std::make_pair(0, "north")),
160 east = interfaces.at(
161 std::make_pair(1, "east")),
162 south = interfaces.at(
163 std::make_pair(2, "south"));
164 BOOST_TEST(north->m_leftInterface == interface1);
165 BOOST_TEST(north->m_rightInterface == interface2);
166 BOOST_TEST(east->m_leftInterface == interface2);
167 BOOST_TEST(east->m_rightInterface == interface3);
168 BOOST_TEST(south->m_leftInterface == interface3);
169 BOOST_TEST(south->m_rightInterface == interface1);
170}
const InterfaceCollection & GetInterfaces() const
Definition: Movement.h:68
void AddInterface(std::string name, InterfaceShPtr left, InterfaceShPtr right)
Add pair of interfaces to this data.
SpatialDomains::InterfaceShPtr CreateInterface(int interfaceID, std::vector< int > compositeIDs)
Produce dummy Interface objects, containing empty domain pointers.
std::shared_ptr< InterfacePair > InterfacePairShPtr
std::shared_ptr< Interface > InterfaceShPtr
std::map< std::pair< int, std::string >, InterfacePairShPtr > InterfaceCollection
Definition: Movement.h:49

References Nektar::SpatialDomains::Movement::AddInterface(), CreateInterface(), and Nektar::SpatialDomains::Movement::GetInterfaces().

◆ BOOST_AUTO_TEST_CASE() [2/3]

Nektar::MovementTests::BOOST_AUTO_TEST_CASE ( TestAddGetZones  )

Definition at line 126 of file TestMovement.cpp.

127{
133 0, interpreter),
134 zone2 = CreateZone(
136 2, interpreter);
137 m.AddZone(zone1);
138 m.AddZone(zone2);
139 std::map<int, SpatialDomains::ZoneBaseShPtr> zones = m.GetZones();
140 BOOST_TEST(zones.size() == 2);
141 BOOST_TEST(zones.at(0).get() == zone1.get());
142 BOOST_TEST(zones.at(2).get() == zone2.get());
143 BOOST_TEST(zones.at(0) == zone1);
144 BOOST_TEST(zones.at(2) == zone2);
145}
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
const std::map< int, ZoneBaseShPtr > & GetZones() const
Definition: Movement.h:73
void AddZone(ZoneBaseShPtr zone)
Add a zone object to this Movement data.
std::shared_ptr< Interpreter > InterpreterSharedPtr
Definition: Interpreter.h:322
SpatialDomains::ZoneBaseShPtr CreateZone(SpatialDomains::MovementType type, int zoneID, int domainID, LibUtilities::InterpreterSharedPtr interpreter)
Produce dummy Zone objects, containing empty domain pointers.
std::shared_ptr< ZoneBase > ZoneBaseShPtr
Definition: Zones.h:170

References Nektar::SpatialDomains::Movement::AddZone(), Nektar::MemoryManager< DataType >::AllocateSharedPtr(), CreateZone(), Nektar::SpatialDomains::eFixed, Nektar::SpatialDomains::eRotate, and Nektar::SpatialDomains::Movement::GetZones().

◆ BOOST_AUTO_TEST_CASE() [3/3]

Nektar::MovementTests::BOOST_AUTO_TEST_CASE ( TestWriteMovement  )

Definition at line 172 of file TestMovement.cpp.

173{
177 m.AddZone(
180 interpreter));
181 m.AddZone(
183 m.AddInterface("north", CreateInterface(0, {0}),
184 CreateInterface(1, {1, 2, 3, 4}));
185 m.AddInterface("south", CreateInterface(2, {5, 6, 7, 8, 9, 11}),
186 CreateInterface(3, {12}));
187
189 std::map<int, SpatialDomains::ZoneBaseShPtr> zones = m.GetZones();
190 TiXmlElement *nektar = new TiXmlElement("NEKTAR");
191 m.WriteMovement(nektar);
192 TiXmlElement *movement = nektar->FirstChildElement("MOVEMENT");
193 BOOST_TEST(movement);
194 TiXmlElement *xmlZones = movement->FirstChildElement("ZONES"),
195 *xmlInterfaces = movement->FirstChildElement("INTERFACES");
196 BOOST_TEST(xmlZones);
197 BOOST_TEST(xmlInterfaces);
198
199 int id, err;
200 std::string attr;
201 std::vector<NekDouble> vec;
202
203 // Check the fixed zone
204 TiXmlElement *zone = xmlZones->FirstChildElement();
205 BOOST_TEST(zone->Value() == "F");
206 err = zone->QueryIntAttribute("ID", &id);
207 BOOST_TEST(err == TIXML_SUCCESS);
208 BOOST_TEST(id == 0);
209 err = zone->QueryStringAttribute("DOMAIN", &attr);
210 BOOST_TEST(err == TIXML_SUCCESS);
211 BOOST_TEST(attr == "D[10]");
212
213 // Check the rotating zone
214 zone = xmlZones->IterateChildren(zone)->ToElement();
215 BOOST_TEST(zone->Value() == "R");
216 err = zone->QueryIntAttribute("ID", &id);
217 BOOST_TEST(err == TIXML_SUCCESS);
218 BOOST_TEST(id == 1);
219 err = zone->QueryStringAttribute("DOMAIN", &attr);
220 BOOST_TEST(err == TIXML_SUCCESS);
221 BOOST_TEST(attr == "D[11]");
222 err = zone->QueryStringAttribute("ORIGIN", &attr);
223 BOOST_TEST(err == TIXML_SUCCESS);
224 ParseUtils::GenerateVector(attr, vec);
225 BOOST_TEST(NekPoint<NekDouble>(vec.at(0), vec.at(1), vec.at(2)) == origin);
226 err = zone->QueryStringAttribute("AXIS", &attr);
227 BOOST_TEST(err == TIXML_SUCCESS);
228 DNekVec xmlAxis(attr);
229 BOOST_TEST(xmlAxis.GetDimension() == 3);
230 for (int i = 0; i < 3; i++)
231 {
232 BOOST_TEST(xmlAxis[i] == axis[i]);
233 }
234 err = zone->QueryStringAttribute("ANGVEL", &attr);
235 BOOST_TEST(err == TIXML_SUCCESS);
236 BOOST_TEST(attr == angVelStr);
237
238 // Check the translating zone
239 zone = xmlZones->IterateChildren(zone)->ToElement();
240 BOOST_TEST(zone->Value() == "T");
241 err = zone->QueryIntAttribute("ID", &id);
242 BOOST_TEST(err == TIXML_SUCCESS);
243 BOOST_TEST(id == 3);
244 err = zone->QueryStringAttribute("DOMAIN", &attr);
245 BOOST_TEST(err == TIXML_SUCCESS);
246 BOOST_TEST(attr == "D[13]");
247 err = zone->QueryStringAttribute("XVELOCITY", &attr);
248 BOOST_TEST(err == TIXML_SUCCESS);
249 BOOST_TEST(attr == velocityStr[0]);
250 err = zone->QueryStringAttribute("YVELOCITY", &attr);
251 BOOST_TEST(err == TIXML_SUCCESS);
252 BOOST_TEST(attr == velocityStr[1]);
253 err = zone->QueryStringAttribute("ZVELOCITY", &attr);
254 BOOST_TEST(err == TIXML_SUCCESS);
255 BOOST_TEST(attr == velocityStr[2]);
256 err = zone->QueryStringAttribute("XDISPLACEMENT", &attr);
257 BOOST_TEST(err == TIXML_SUCCESS);
258 BOOST_TEST(attr == displacementStr[0]);
259 err = zone->QueryStringAttribute("YDISPLACEMENT", &attr);
260 BOOST_TEST(err == TIXML_SUCCESS);
261 BOOST_TEST(attr == displacementStr[1]);
262 err = zone->QueryStringAttribute("ZDISPLACEMENT", &attr);
263 BOOST_TEST(err == TIXML_SUCCESS);
264 BOOST_TEST(attr == displacementStr[2]);
265
266 BOOST_TEST(xmlZones->LastChild() == zone);
267
268 TiXmlElement *intr;
269
270 // Check the north interface
271 TiXmlElement *interface = xmlInterfaces->FirstChildElement();
272 BOOST_TEST(interface->Value() == "INTERFACE");
273 err = interface->QueryStringAttribute("NAME", &attr);
274 BOOST_TEST(err == TIXML_SUCCESS);
275 BOOST_TEST(attr == "north");
276 intr = interface->FirstChildElement();
277 BOOST_TEST(intr->Value() == "L");
278 err = intr->QueryIntAttribute("ID", &id);
279 BOOST_TEST(err == TIXML_SUCCESS);
280 BOOST_TEST(id == 0);
281 err = intr->QueryStringAttribute("BOUNDARY", &attr);
282 BOOST_TEST(err == TIXML_SUCCESS);
283 BOOST_TEST(attr == "C[0]");
284 intr = interface->IterateChildren(intr)->ToElement();
285 BOOST_TEST(intr->Value() == "R");
286 err = intr->QueryIntAttribute("ID", &id);
287 BOOST_TEST(err == TIXML_SUCCESS);
288 BOOST_TEST(id == 1);
289 err = intr->QueryStringAttribute("BOUNDARY", &attr);
290 BOOST_TEST(err == TIXML_SUCCESS);
291 BOOST_TEST(attr == "C[1-4]");
292
293 // Check the south interface
294 interface = xmlInterfaces->IterateChildren(interface)->ToElement();
295 BOOST_TEST(interface->Value() == "INTERFACE");
296 err = interface->QueryStringAttribute("NAME", &attr);
297 BOOST_TEST(err == TIXML_SUCCESS);
298 BOOST_TEST(attr == "south");
299 intr = interface->FirstChildElement();
300 BOOST_TEST(intr->Value() == "L");
301 err = intr->QueryIntAttribute("ID", &id);
302 BOOST_TEST(err == TIXML_SUCCESS);
303 BOOST_TEST(id == 2);
304 err = intr->QueryStringAttribute("BOUNDARY", &attr);
305 BOOST_TEST(err == TIXML_SUCCESS);
306 BOOST_TEST(attr == "C[5-9,11]");
307 intr = interface->IterateChildren(intr)->ToElement();
308 BOOST_TEST(intr->Value() == "R");
309 err = intr->QueryIntAttribute("ID", &id);
310 BOOST_TEST(err == TIXML_SUCCESS);
311 BOOST_TEST(id == 3);
312 err = intr->QueryStringAttribute("BOUNDARY", &attr);
313 BOOST_TEST(err == TIXML_SUCCESS);
314 BOOST_TEST(attr == "C[12]");
315
316 BOOST_TEST(xmlInterfaces->LastChild() == interface);
317
318 delete nektar;
319}
void WriteMovement(TiXmlElement *root)
Write the MOVEMENT section of the XML file.
std::vector< std::string > velocityStr
const NekPoint< NekDouble > origin
const std::string angVelStr
std::vector< std::string > displacementStr
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48

References Nektar::SpatialDomains::Movement::AddInterface(), Nektar::SpatialDomains::Movement::AddZone(), Nektar::MemoryManager< DataType >::AllocateSharedPtr(), angVelStr, axis, CreateInterface(), CreateZone(), displacementStr, Nektar::SpatialDomains::eFixed, Nektar::SpatialDomains::eRotate, Nektar::SpatialDomains::eTranslate, Nektar::ParseUtils::GenerateVector(), Nektar::NekVector< DataType >::GetDimension(), Nektar::SpatialDomains::Movement::GetInterfaces(), Nektar::SpatialDomains::Movement::GetZones(), origin, velocityStr, and Nektar::SpatialDomains::Movement::WriteMovement().

◆ CreateInterface()

SpatialDomains::InterfaceShPtr Nektar::MovementTests::CreateInterface ( int  interfaceID,
std::vector< int >  compositeIDs 
)

Produce dummy Interface objects, containing empty domain pointers.

Definition at line 112 of file TestMovement.cpp.

114{
116 for (auto &id : compositeIDs)
117 {
118 edge[id] =
120 }
123 interfaceID, edge, false));
124}
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:136

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

Referenced by BOOST_AUTO_TEST_CASE().

◆ CreateSession()

LibUtilities::SessionReaderSharedPtr Nektar::MovementTests::CreateSession ( )

Definition at line 46 of file TestMovement.cpp.

47{
48 char arg = ' ';
49 char *args[] = {&arg};
50 return LibUtilities::SessionReader::CreateInstance(1, args);
51}

References Nektar::LibUtilities::SessionReader::CreateInstance().

◆ CreateZone()

SpatialDomains::ZoneBaseShPtr Nektar::MovementTests::CreateZone ( SpatialDomains::MovementType  type,
int  zoneID,
int  domainID,
LibUtilities::InterpreterSharedPtr  interpreter 
)

Produce dummy Zone objects, containing empty domain pointers.

Definition at line 61 of file TestMovement.cpp.

64{
66 switch (type)
67 {
69 {
72 zoneID, domainID, domain, 3));
73 }
74 break;
76 {
78 std::make_shared<LibUtilities::Equation>(interpreter,
79 angVelStr);
82 zoneID, domainID, domain, 3, origin, axis, angularVelEqn));
83 }
84 break;
86 {
87 Array<OneD, LibUtilities::EquationSharedPtr> velocityEqns(3);
88 Array<OneD, LibUtilities::EquationSharedPtr> displacementEqns(3);
89 for (int i = 0; i < 3; ++i)
90 {
91 velocityEqns[i] = std::make_shared<LibUtilities::Equation>(
92 interpreter, velocityStr[i]);
93 displacementEqns[i] = std::make_shared<LibUtilities::Equation>(
94 interpreter, displacementStr[i]);
95 }
96
98 MemoryManager<SpatialDomains::ZoneTranslate>::AllocateSharedPtr(
99 zoneID, domainID, domain, 3, velocityEqns,
100 displacementEqns));
101 }
102 break;
103 default:
104 {
105 }
106 break;
107 }
108 return nullptr;
109}
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
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

References angVelStr, axis, displacementStr, Nektar::SpatialDomains::eFixed, Nektar::SpatialDomains::eRotate, Nektar::SpatialDomains::eTranslate, origin, and velocityStr.

Referenced by BOOST_AUTO_TEST_CASE().

Variable Documentation

◆ angVelStr

const std::string Nektar::MovementTests::angVelStr = "0.1*t"

Definition at line 53 of file TestMovement.cpp.

Referenced by BOOST_AUTO_TEST_CASE(), and CreateZone().

◆ axis

const DNekVec Nektar::MovementTests::axis = {1., 2., 3.}

◆ displacementStr

std::vector<std::string> Nektar::MovementTests::displacementStr = {"1.0", "2.0", "3.0"}

◆ origin

const NekPoint<NekDouble> Nektar::MovementTests::origin = {1., 2., 3.}

◆ velocityStr

std::vector<std::string> Nektar::MovementTests::velocityStr = {"1.0", "2.0", "3.0"}

◆ xEqnStr

const std::string Nektar::MovementTests::xEqnStr = "0.1*x - 0.1*t"

Definition at line 53 of file TestMovement.cpp.

◆ yEqnStr

const std::string Nektar::MovementTests::yEqnStr = "0.1*y^2 - x"

Definition at line 54 of file TestMovement.cpp.

◆ zEqnStr

const std::string Nektar::MovementTests::zEqnStr = "sqrt(t)"

Definition at line 54 of file TestMovement.cpp.