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)"
 
const NekPoint< NekDoubleorigin = {1., 2., 3.}
 
const DNekVec axis = {1., 2., 3.}
 
const std::vector< NekDoublevelocity = {1., 1., 2.}
 

Function Documentation

◆ BOOST_AUTO_TEST_CASE() [1/3]

Nektar::MovementTests::BOOST_AUTO_TEST_CASE ( TestAddGetInterfaces  )

Definition at line 149 of file TestMovement.cpp.

150{
153 interface2 = CreateInterface(1, {1, 2, 3}),
154 interface3 = CreateInterface(3, {4});
155 m.AddInterface("north", interface1, interface2);
156 m.AddInterface("east", interface2, interface3);
157 m.AddInterface("south", interface3, interface1);
159 BOOST_TEST(interfaces.size() == 3);
160 SpatialDomains::InterfacePairShPtr north = interfaces.at(
161 std::make_pair(0, "north")),
162 east = interfaces.at(
163 std::make_pair(1, "east")),
164 south = interfaces.at(
165 std::make_pair(2, "south"));
166 BOOST_TEST(north->m_leftInterface == interface1);
167 BOOST_TEST(north->m_rightInterface == interface2);
168 BOOST_TEST(east->m_leftInterface == interface2);
169 BOOST_TEST(east->m_rightInterface == interface3);
170 BOOST_TEST(south->m_leftInterface == interface3);
171 BOOST_TEST(south->m_rightInterface == interface1);
172}
const InterfaceCollection & GetInterfaces() const
Definition: Movement.h:69
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 128 of file TestMovement.cpp.

129{
135 0, interpreter),
136 zone2 = CreateZone(
138 2, interpreter);
139 m.AddZone(zone1);
140 m.AddZone(zone2);
141 std::map<int, SpatialDomains::ZoneBaseShPtr> zones = m.GetZones();
142 BOOST_TEST(zones.size() == 2);
143 BOOST_TEST(zones.at(0).get() == zone1.get());
144 BOOST_TEST(zones.at(2).get() == zone2.get());
145 BOOST_TEST(zones.at(0) == zone1);
146 BOOST_TEST(zones.at(2) == zone2);
147}
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:74
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:144

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 174 of file TestMovement.cpp.

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

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

◆ CreateInterface()

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

Produce dummy Interface objects, containing empty domain pointers.

Definition at line 114 of file TestMovement.cpp.

116{
118 for (auto &id : compositeIDs)
119 {
120 edge[id] =
122 }
125 edge));
126}
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 60 of file TestMovement.cpp.

63{
65 switch (type)
66 {
68 {
71 zoneID, domainID, domain, 3));
72 }
73 break;
75 {
77 std::make_shared<LibUtilities::Equation>(interpreter,
78 angVelStr);
81 zoneID, domainID, domain, 3, origin, axis, angularVelEqn));
82 }
83 break;
85 {
87 MemoryManager<SpatialDomains::ZoneTranslate>::AllocateSharedPtr(
88 zoneID, domainID, domain, 3, velocity));
89 }
90 break;
92 {
94 xEqn = std::make_shared<LibUtilities::Equation>(interpreter,
95 xEqnStr),
96 yEqn = std::make_shared<LibUtilities::Equation>(interpreter,
97 yEqnStr),
98 zEqn = std::make_shared<LibUtilities::Equation>(interpreter,
99 zEqnStr);
101 MemoryManager<SpatialDomains::ZonePrescribe>::AllocateSharedPtr(
102 zoneID, domainID, domain, 3, xEqn, yEqn, zEqn));
103 }
104 break;
105 default:
106 {
107 }
108 break;
109 }
110 return nullptr;
111}
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
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

References angVelStr, axis, Nektar::SpatialDomains::eFixed, Nektar::SpatialDomains::ePrescribe, Nektar::SpatialDomains::eRotate, Nektar::SpatialDomains::eTranslate, origin, velocity, xEqnStr, yEqnStr, and zEqnStr.

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

◆ origin

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

◆ velocity

const std::vector<NekDouble> Nektar::MovementTests::velocity = {1., 1., 2.}

Definition at line 57 of file TestMovement.cpp.

Referenced by BOOST_AUTO_TEST_CASE(), Nektar::SolverUtils::FilterAeroForces::CalculateForces(), CreateZone(), Nektar::IncNavierStokes::EvaluateAdvectionTerms(), Nektar::MMFAdvection::EvaluateAdvectionVelocity(), Nektar::VariableConverter::GetAbsoluteVelocity(), Nektar::NonlinearPeregrine::GetFluxVector(), Nektar::NonlinearSWE::GetFluxVector(), Nektar::CompressibleFlowSystem::GetFluxVector(), Nektar::CompressibleFlowSystem::GetFluxVectorDeAlias(), Nektar::SolverUtils::FluidInterface::GetVelocity(), Nektar::VariableConverter::GetVelocityVector(), Nektar::LinearSWE::GetVelocityVector(), Nektar::NonlinearPeregrine::GetVelocityVector(), Nektar::NonlinearSWE::GetVelocityVector(), Nektar::MultiRegions::ExpList::LinearAdvectionReactionSolve(), main(), OUTPUT(), Nektar::SpatialDomains::Movement::ReadZones(), Nektar::AdjointAdvection::v_Advect(), Nektar::AlternateSkewAdvection::v_Advect(), Nektar::LinearisedAdvection::v_Advect(), Nektar::NavierStokesAdvection::v_Advect(), Nektar::SkewSymmetricAdvection::v_Advect(), Nektar::CompressibleFlowSystem::v_ExtraFldOutput(), Nektar::NavierStokesCFE::v_ExtraFldOutput(), Nektar::AcousticSystem::v_GetMaxStdVelocity(), Nektar::CompressibleFlowSystem::v_GetMaxStdVelocity(), Nektar::CompressibleFlowSystem::v_GetVelocity(), Nektar::IncNavierStokes::v_GetVelocity(), Nektar::SolverUtils::FileSolution::v_GetVelocity(), Nektar::CoupledLinearNS::v_InitObject(), Nektar::MultiRegions::ContField::v_LinearAdvectionReactionSolve(), Nektar::FieldUtils::ProcessWSS::v_Process(), Nektar::VCSMapping::v_SetUpPressureForcing(), Nektar::VCSImplicit::v_SetUpViscousForcing(), and ZoneTranslate_Init().

◆ xEqnStr

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

Definition at line 53 of file TestMovement.cpp.

Referenced by BOOST_AUTO_TEST_CASE(), and CreateZone().

◆ yEqnStr

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

Definition at line 54 of file TestMovement.cpp.

Referenced by BOOST_AUTO_TEST_CASE(), and CreateZone().

◆ zEqnStr

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

Definition at line 54 of file TestMovement.cpp.

Referenced by BOOST_AUTO_TEST_CASE(), and CreateZone().