PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
singlebody_test.cpp
Go to the documentation of this file.
1
8#include <vector>
9
10#include <petsc.h>
11
12#include <gtest/gtest.h>
13#include <yaml-cpp/yaml.h>
14
15#include <petibm/mesh.h>
16#include <petibm/singlebody.h>
17
18using namespace petibm;
19
20class SingleBodyTest : public ::testing::Test
21{
22protected:
23 SingleBodyTest(){};
24
25 virtual ~SingleBodyTest(){};
26
27 virtual void SetUp()
28 {
29 using namespace YAML;
30 YAML::Node config;
31
32 config["mesh"].push_back(Node(NodeType::Map));
33 config["mesh"][0]["direction"] = "x";
34 config["mesh"][1]["direction"] = "y";
35 for (unsigned int i = 0; i < 2; ++i)
36 {
37 config["mesh"][i]["start"] = "0.1";
38 config["mesh"][i]["subDomains"].push_back(Node(NodeType::Map));
39 config["mesh"][i]["subDomains"][0]["end"] = 1.0;
40 config["mesh"][i]["subDomains"][0]["cells"] = 10;
41 config["mesh"][i]["subDomains"][0]["stretchRatio"] = 1.0;
42 }
43
44 config["flow"] = YAML::Node(NodeType::Map);
45 config["flow"]["boundaryConditions"].push_back(Node(NodeType::Map));
46 config["flow"]["boundaryConditions"][0]["location"] = "xMinus";
47 config["flow"]["boundaryConditions"][1]["location"] = "xPlus";
48 config["flow"]["boundaryConditions"][2]["location"] = "yMinus";
49 config["flow"]["boundaryConditions"][3]["location"] = "yPlus";
50
51 for (unsigned int i = 0; i < 4; ++i)
52 {
53 config["flow"]["boundaryConditions"][i]["u"][0] = "DIRICHLET";
54 config["flow"]["boundaryConditions"][i]["u"][1] = 0.0;
55 config["flow"]["boundaryConditions"][i]["v"][0] = "DIRICHLET";
56 config["flow"]["boundaryConditions"][i]["v"][1] = 0.0;
57 }
58 config["flow"]["boundaryConditions"][3]["u"][1] = 1.0;
59
60 config["directory"] = "body/";
61 config["bodies"][0]["type"] = "points";
62 config["bodies"][0]["file"] = "body2d.txt";
63
64 // create 2D mesh and body
65 petibm::mesh::createMesh(PETSC_COMM_WORLD, config, mesh2d);
66 petibm::body::createSingleBody(PETSC_COMM_WORLD, 2, "points",
67 "body2d01", "body/body2d.txt", body2d);
68
69 // create 3D mesh and body
70 config["mesh"][2]["direction"] = "z";
71 config["mesh"][2]["start"] = "0.1";
72 config["mesh"][2]["subDomains"].push_back(Node(NodeType::Map));
73 config["mesh"][2]["subDomains"][0]["end"] = 1.0;
74 config["mesh"][2]["subDomains"][0]["cells"] = 10;
75 config["mesh"][2]["subDomains"][0]["stretchRatio"] = 1.0;
76 for (unsigned int i = 0; i < 4; ++i)
77 {
78 config["flow"]["boundaryConditions"][i]["w"][0] = "DIRICHLET";
79 config["flow"]["boundaryConditions"][i]["w"][1] = 0.0;
80 }
81 config["flow"]["boundaryConditions"][4]["location"] = "zMinus";
82 config["flow"]["boundaryConditions"][5]["location"] = "zPlus";
83
84 for (unsigned int i = 4; i < 6; ++i)
85 {
86 config["flow"]["boundaryConditions"][i]["u"][0] = "PERIODIC";
87 config["flow"]["boundaryConditions"][i]["u"][1] = 0.0;
88 config["flow"]["boundaryConditions"][i]["v"][0] = "PERIODIC";
89 config["flow"]["boundaryConditions"][i]["v"][1] = 0.0;
90 config["flow"]["boundaryConditions"][i]["w"][0] = "PERIODIC";
91 config["flow"]["boundaryConditions"][i]["w"][1] = 0.0;
92 }
93 config["bodies"][0]["file"] = "body3d.txt";
94
95 petibm::mesh::createMesh(PETSC_COMM_WORLD, config, mesh3d);
96 petibm::body::createSingleBody(PETSC_COMM_WORLD, 3, "points",
97 "body3d01", "body/body3d.txt", body3d);
98 };
99
100 virtual void TearDown(){};
101
102 type::SingleBody body2d, body3d;
103 type::Mesh mesh2d, mesh3d;
104
105}; // SingleBodyTest
106
107TEST_F(SingleBodyTest, initWithFilePath2D)
108{
109 ASSERT_EQ(2, body2d->dim);
110 ASSERT_EQ("body2d01", body2d->name);
111 ASSERT_EQ(4, body2d->nPts);
112 std::vector<PetscReal> xCoords = {0.25, 0.75, 0.75, 0.25},
113 yCoords = {0.25, 0.25, 0.75, 0.75};
114 for (unsigned int i = 0; i < 4; i++)
115 {
116 ASSERT_EQ(xCoords[i], body2d->coords[i][0]);
117 ASSERT_EQ(yCoords[i], body2d->coords[i][1]);
118 }
119 PetscMPIInt size;
120 MPI_Comm_size(PETSC_COMM_WORLD, &size);
121 if (size == 1) ASSERT_EQ(body2d->nPts, body2d->nLclPts);
122}
123
124TEST_F(SingleBodyTest, initWithFilePath3D)
125{
126 ASSERT_EQ(3, body3d->dim);
127 ASSERT_EQ("body3d01", body3d->name);
128 ASSERT_EQ(8, body3d->nPts);
129 std::vector<PetscReal> xCoords = {0.25, 0.75, 0.75, 0.25,
130 0.25, 0.75, 0.75, 0.25},
131 yCoords = {0.25, 0.25, 0.75, 0.75,
132 0.25, 0.25, 0.75, 0.75},
133 zCoords = {0.25, 0.25, 0.25, 0.25,
134 0.75, 0.75, 0.75, 0.75};
135 for (unsigned int i = 0; i < 4; i++)
136 {
137 ASSERT_EQ(xCoords[i], body3d->coords[i][0]);
138 ASSERT_EQ(yCoords[i], body3d->coords[i][1]);
139 ASSERT_EQ(zCoords[i], body3d->coords[i][2]);
140 }
141 PetscMPIInt size;
142 MPI_Comm_size(PETSC_COMM_WORLD, &size);
143 if (size == 1) ASSERT_EQ(body3d->nPts, body3d->nLclPts);
144}
145
146TEST_F(SingleBodyTest, findProc2D)
147{
148 PetscMPIInt size, index;
149 MPI_Comm_size(PETSC_COMM_WORLD, &size);
150 if (size == 1)
151 for (int i = 0; i < body2d->nPts; i++)
152 {
153 body2d->findProc(i, index);
154 ASSERT_EQ(0, index);
155 }
156}
157
158TEST_F(SingleBodyTest, findProc3D)
159{
160 PetscMPIInt size, index;
161 MPI_Comm_size(PETSC_COMM_WORLD, &size);
162 if (size == 1)
163 for (int i = 0; i < body3d->nPts; i++)
164 {
165 body3d->findProc(i, index);
166 ASSERT_EQ(0, index);
167 }
168}
169
170TEST_F(SingleBodyTest, getGlobalIndex2D)
171{
172 PetscInt globalIndex, counter = 0;
173 PetscInt ndof = 2;
174 for (int i = 0; i < body2d->nPts; i++)
175 for (int d = 0; d < ndof; d++)
176 {
177 body2d->getGlobalIndex(i, d, globalIndex);
178 ASSERT_EQ(counter, globalIndex);
179 MatStencil stencil = {0, 0, i, d};
180 body2d->getGlobalIndex(stencil, globalIndex);
181 ASSERT_EQ(counter, globalIndex);
182 counter++;
183 }
184}
185
186TEST_F(SingleBodyTest, getGlobalIndex3D)
187{
188 PetscInt globalIndex, counter = 0;
189 PetscInt ndof = 3;
190 for (int i = 0; i < body3d->nPts; i++)
191 for (int d = 0; d < ndof; d++)
192 {
193 body3d->getGlobalIndex(i, d, globalIndex);
194 ASSERT_EQ(counter, globalIndex);
195 MatStencil stencil = {0, 0, i, d};
196 body3d->getGlobalIndex(stencil, globalIndex);
197 ASSERT_EQ(counter, globalIndex);
198 counter++;
199 }
200}
201
202TEST_F(SingleBodyTest, calculateAvgForces2D)
203{
204 Vec f;
205 DMCreateGlobalVector(body2d->da, &f);
206 VecSet(f, 1.0);
207 type::RealVec1D avg(2);
208 body2d->calculateAvgForces(f, avg);
209 ASSERT_EQ(body2d->dim, (int)avg.size());
210 for (unsigned int d = 0; d < avg.size(); d++)
211 ASSERT_EQ(body2d->nPts * -1.0, avg[d]);
212 VecDestroy(&f);
213}
214
215TEST_F(SingleBodyTest, calculateAvgForces3D)
216{
217 Vec f;
218 DMCreateGlobalVector(body3d->da, &f);
219 VecSet(f, 1.0);
220 type::RealVec1D avg(3);
221 body3d->calculateAvgForces(f, avg);
222 ASSERT_EQ(body3d->dim, (int)avg.size());
223 for (unsigned int d = 0; d < avg.size(); d++)
224 ASSERT_EQ(body3d->nPts * -1.0, avg[d]);
225 VecDestroy(&f);
226}
227
228// Run all tests
229int main(int argc, char **argv)
230{
231 PetscErrorCode ierr, status;
232
233 ::testing::InitGoogleTest(&argc, argv);
234 ierr = PetscInitialize(&argc, &argv, nullptr, nullptr); CHKERRQ(ierr);
235 status = RUN_ALL_TESTS();
236 ierr = PetscFinalize(); CHKERRQ(ierr);
237
238 return status;
239} // main
PetscErrorCode createSingleBody(const MPI_Comm &comm, const PetscInt &dim, const std::string &type, const std::string &name, const std::string &filePath, type::SingleBody &body)
Factory function to create a single body.
Definition: singlebody.cpp:83
std::shared_ptr< body::SingleBodyBase > SingleBody
Definition of type::SingleBody.
Definition: singlebody.h:206
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
Definition: mesh.h:348
PetscErrorCode createMesh(const MPI_Comm &comm, const YAML::Node &node, type::Mesh &mesh)
Factory function for creating a Mesh object.
Definition: mesh.cpp:86
std::vector< PetscReal > RealVec1D
1D std::vector holding PetscReal.
Definition: type.h:140
Prototype of mesh::MeshBase, type::Mesh, and factory function.
Definition: parser.cpp:45
A toolbox for building flow solvers.
Definition: bodypack.h:52
body::SingleBodyBase, type::SingleBody factory function.
TEST_F(SingleBodyTest, initWithFilePath2D)
int main(int argc, char **argv)