PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
cartesianmesh2d_yperiodic.cpp
Go to the documentation of this file.
1
8#include <petsc.h>
9
10#include <gtest/gtest.h>
11#include <yaml-cpp/yaml.h>
12
13#include <petibm/mesh.h>
14
15class CartesainMeshTest2D_YPeriodic : public ::testing::Test
16{
17protected:
18 CartesainMeshTest2D_YPeriodic(){};
19
20 virtual ~CartesainMeshTest2D_YPeriodic(){};
21
22 static void SetUpTestCase()
23 {
24 using namespace YAML;
25
26 Node config;
27
28 config["mesh"].push_back(Node(NodeType::Map));
29
30 config["mesh"][0]["direction"] = "x";
31 config["mesh"][0]["start"] = "0.1";
32 config["mesh"][0]["subDomains"].push_back(Node(NodeType::Map));
33 config["mesh"][0]["subDomains"][0]["end"] = 1.6;
34 config["mesh"][0]["subDomains"][0]["cells"] = 4;
35 config["mesh"][0]["subDomains"][0]["stretchRatio"] = 0.5;
36 config["mesh"][0]["subDomains"][1]["end"] = 1.9;
37 config["mesh"][0]["subDomains"][1]["cells"] = 3;
38 config["mesh"][0]["subDomains"][1]["stretchRatio"] = 1;
39 config["mesh"][0]["subDomains"][2]["end"] = 5.0;
40 config["mesh"][0]["subDomains"][2]["cells"] = 5;
41 config["mesh"][0]["subDomains"][2]["stretchRatio"] = 2.0;
42
43 config["mesh"][1]["direction"] = "y";
44 config["mesh"][1]["start"] = "0.05";
45 config["mesh"][1]["subDomains"].push_back(Node(NodeType::Map));
46 config["mesh"][1]["subDomains"][0]["end"] = 0.8625;
47 config["mesh"][1]["subDomains"][0]["cells"] = 4;
48 config["mesh"][1]["subDomains"][0]["stretchRatio"] = 0.666666666666666;
49 config["mesh"][1]["subDomains"][1]["end"] = 1.1625;
50 config["mesh"][1]["subDomains"][1]["cells"] = 3;
51 config["mesh"][1]["subDomains"][1]["stretchRatio"] = 1.0;
52 config["mesh"][1]["subDomains"][2]["end"] = 1.975;
53 config["mesh"][1]["subDomains"][2]["cells"] = 4;
54 config["mesh"][1]["subDomains"][2]["stretchRatio"] = 1.5;
55
56 config["flow"] = YAML::Node(NodeType::Map);
57 config["flow"]["boundaryConditions"].push_back(Node(NodeType::Map));
58 config["flow"]["boundaryConditions"][0]["location"] = "xMinus";
59 config["flow"]["boundaryConditions"][1]["location"] = "xPlus";
60 config["flow"]["boundaryConditions"][2]["location"] = "yMinus";
61 config["flow"]["boundaryConditions"][3]["location"] = "yPlus";
62
63 for (unsigned int i = 0; i < 2; ++i)
64 {
65 config["flow"]["boundaryConditions"][i]["u"][0] = "DIRICHLET";
66 config["flow"]["boundaryConditions"][i]["u"][1] = 0.0;
67 config["flow"]["boundaryConditions"][i]["v"][0] = "DIRICHLET";
68 config["flow"]["boundaryConditions"][i]["v"][1] = 0.0;
69 }
70
71 for (unsigned int i = 2; i < 4; ++i)
72 {
73 config["flow"]["boundaryConditions"][i]["u"][0] = "PERIODIC";
74 config["flow"]["boundaryConditions"][i]["u"][1] = 0.0;
75 config["flow"]["boundaryConditions"][i]["v"][0] = "PERIODIC";
76 config["flow"]["boundaryConditions"][i]["v"][1] = 0.0;
77 }
78
79 config["flow"]["boundaryConditions"][3]["u"][1] = 1.0;
80
81 PetscErrorCode ierr;
82 PetscInt nProc[3];
83
84 ierr = petibm::mesh::createMesh(PETSC_COMM_WORLD, config, mesh);
85 ASSERT_FALSE(ierr);
86
87 ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_GHOSTED,
88 DM_BOUNDARY_GHOSTED, DMDA_STENCIL_BOX, 12, 11,
89 PETSC_DECIDE, PETSC_DECIDE, 1, 1, nullptr, nullptr,
90 &da[3]);
91 ASSERT_FALSE(ierr);
92 ierr = DMSetUp(da[3]);
93 ASSERT_FALSE(ierr);
94
95 ierr = DMDAGetInfo(da[3], nullptr, nullptr, nullptr, nullptr, &nProc[0],
96 &nProc[1], &nProc[2], nullptr, nullptr, nullptr,
97 nullptr, nullptr, nullptr);
98 ASSERT_FALSE(ierr);
99
100 ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_GHOSTED,
101 DM_BOUNDARY_GHOSTED, DMDA_STENCIL_BOX, 11, 11,
102 nProc[0], nProc[1], 1, 1, nullptr, nullptr, &da[0]);
103 ASSERT_FALSE(ierr);
104 ierr = DMSetUp(da[0]);
105 ASSERT_FALSE(ierr);
106
107 ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_GHOSTED,
108 DM_BOUNDARY_GHOSTED, DMDA_STENCIL_BOX, 12, 11,
109 nProc[0], nProc[1], 1, 1, nullptr, nullptr, &da[1]);
110 ASSERT_FALSE(ierr);
111 ierr = DMSetUp(da[1]);
112 ASSERT_FALSE(ierr);
113 };
114
115 virtual void SetUp(){};
116
117 virtual void TearDown(){};
118
119 static void TearDownTestCase()
120 {
121 PetscErrorCode ierr;
122 ierr = DMDestroy(&da[0]);
123 ASSERT_FALSE(ierr);
124 ierr = DMDestroy(&da[1]);
125 ASSERT_FALSE(ierr);
126 ierr = DMDestroy(&da[3]);
127 ASSERT_FALSE(ierr);
128 };
129
130 static petibm::type::Mesh mesh;
131 static DM da[5];
132}; // CartesianMeshTest
133
134petibm::type::Mesh CartesainMeshTest2D_YPeriodic::mesh = nullptr;
135DM CartesainMeshTest2D_YPeriodic::da[5]{nullptr, nullptr, nullptr, nullptr,
136 nullptr};
137
138// test dimension
139TEST_F(CartesainMeshTest2D_YPeriodic, check_dim) { ASSERT_EQ(2, mesh->dim); }
140
141// check min
142TEST_F(CartesainMeshTest2D_YPeriodic, check_min)
143{
144 ASSERT_DOUBLE_EQ(0.1, mesh->min[0]);
145 ASSERT_DOUBLE_EQ(0.05, mesh->min[1]);
146 ASSERT_DOUBLE_EQ(0.0, mesh->min[2]);
147}
148
149// check max
150TEST_F(CartesainMeshTest2D_YPeriodic, check_max)
151{
152 ASSERT_DOUBLE_EQ(5.0, mesh->max[0]);
153 ASSERT_DOUBLE_EQ(1.975, mesh->max[1]);
154 ASSERT_DOUBLE_EQ(1.0, mesh->max[2]);
155}
156
157// check n
158TEST_F(CartesainMeshTest2D_YPeriodic, check_n)
159{
161 {11, 11, 1}, {12, 11, 1}, {1, 1, 1}, {12, 11, 1}, {13, 12, 1}};
162
163 for (unsigned int f = 0; f < 5; ++f)
164 for (unsigned int d = 0; d < 3; ++d) ASSERT_EQ(n[f][d], mesh->n[f][d]);
165}
166
167// check periodic
168TEST_F(CartesainMeshTest2D_YPeriodic, check_periodic)
169{
170 petibm::type::BoolVec2D periodic{{PETSC_FALSE, PETSC_TRUE, PETSC_FALSE},
171 {PETSC_FALSE, PETSC_TRUE, PETSC_FALSE},
172 {PETSC_FALSE, PETSC_FALSE, PETSC_FALSE}};
173
174 for (unsigned int f = 0; f < 3; ++f)
175 for (unsigned int d = 0; d < 3; ++d)
176 ASSERT_EQ(periodic[f][d], mesh->periodic[f][d]);
177}
178
179// check coord
180TEST_F(CartesainMeshTest2D_YPeriodic, check_coord)
181{
182 petibm::type::RealVec3D coordTrue{
183 {{0.1, 0.9, 1.3, 1.5, 1.6, 1.7, 1.8, 1.9, 2.0, 2.2, 2.6, 3.4, 5.0},
184 {-0.11875, 0.21875, 0.5, 0.6875, 0.8125, 0.9125, 1.0125, 1.1125,
185 1.2125, 1.3375, 1.525, 1.80625, 2.14375},
186 {0.0}},
187 {{-0.3, 0.5, 1.1, 1.4, 1.55, 1.65, 1.75, 1.85, 1.95, 2.1, 2.4, 3.0, 4.2,
188 5.8},
189 {0.05, 0.3875, 0.6125, 0.7625, 0.8625, 0.9625, 1.0625, 1.1625, 1.2625,
190 1.4125, 1.6375, 1.975, 2.3125},
191 {0.0}},
192 {{0.0}, {0.0}, {0.0}},
193 {{0.5, 1.1, 1.4, 1.55, 1.65, 1.75, 1.85, 1.95, 2.1, 2.4, 3.0, 4.2},
194 {0.21875, 0.5, 0.6875, 0.8125, 0.9125, 1.0125, 1.1125, 1.2125, 1.3375,
195 1.525, 1.80625},
196 {0.0}},
197 {{0.1, 0.9, 1.3, 1.5, 1.6, 1.7, 1.8, 1.9, 2.0, 2.2, 2.6, 3.4, 5.0},
198 {0.05, 0.3875, 0.6125, 0.7625, 0.8625, 0.9625, 1.0625, 1.1625, 1.2625,
199 1.4125, 1.6375, 1.975},
200 {0.0}}};
201
203 {&coordTrue[0][0][1], &coordTrue[0][1][1], &coordTrue[0][2][0]},
204 {&coordTrue[1][0][1], &coordTrue[1][1][1], &coordTrue[1][2][0]},
205 {&coordTrue[2][0][0], &coordTrue[2][1][0], &coordTrue[2][2][0]},
206 {&coordTrue[3][0][0], &coordTrue[3][1][0], &coordTrue[3][2][0]},
207 {&coordTrue[4][0][0], &coordTrue[4][1][0], &coordTrue[4][2][0]},
208 };
209
210 for (unsigned int f = 0; f < 2; ++f)
211 {
212 for (unsigned int d = 0; d < 2; ++d)
213 for (int i = -1; i < mesh->n[f][d] + 1; ++i)
214 ASSERT_NEAR(coord[f][d][i], mesh->coord[f][d][i], 1e-12)
215 << "dL incorrect at field " << f << ", direction " << d
216 << ", index " << i << std::endl;
217
218 ASSERT_NEAR(coord[f][2][0], mesh->coord[f][2][0], 1e-12)
219 << "dL incorrect at field " << f << ", direction " << 2
220 << ", index " << 0 << std::endl;
221 }
222
223 for (unsigned int f = 2; f < 5; ++f)
224 for (unsigned int d = 0; d < 3; ++d)
225 for (int i = 0; i < mesh->n[f][d]; ++i)
226 ASSERT_NEAR(coord[f][d][i], mesh->coord[f][d][i], 1e-12)
227 << "dL incorrect at field " << f << ", direction " << d
228 << ", index " << i << std::endl;
229}
230
231// check dL
232TEST_F(CartesainMeshTest2D_YPeriodic, check_dL)
233{
235 {{0.8, 0.6, 0.3, 0.15, 0.1, 0.1, 0.1, 0.1, 0.15, 0.3, 0.6, 1.2, 1.6},
236 {0.3375, 0.3375, 0.225, 0.15, 0.1, 0.1, 0.1, 0.1, 0.1, 0.15, 0.225,
237 0.3375, 0.3375},
238 {1.0}},
239 {{0.8, 0.8, 0.4, 0.2, 0.1, 0.1, 0.1, 0.1, 0.1, 0.2, 0.4, 0.8, 1.6, 1.6},
240 {0.3375, 0.28125, 0.1875, 0.125, 0.1, 0.1, 0.1, 0.1, 0.125, 0.1875,
241 0.28125, 0.3375, 0.28125},
242 {1.0}},
243 {{1.0}, {1.0}, {1.0}},
244 {{0.8, 0.4, 0.2, 0.1, 0.1, 0.1, 0.1, 0.1, 0.2, 0.4, 0.8, 1.6},
245 {0.3375, 0.225, 0.15, 0.1, 0.1, 0.1, 0.1, 0.1, 0.15, 0.225, 0.3375},
246 {1.0}},
247 {{1.0}, {1.0}, {1.0}}};
248
250 {&dLTrue[0][0][1], &dLTrue[0][1][1], &dLTrue[0][2][0]},
251 {&dLTrue[1][0][1], &dLTrue[1][1][1], &dLTrue[1][2][0]},
252 {&dLTrue[2][0][0], &dLTrue[2][1][0], &dLTrue[2][2][0]},
253 {&dLTrue[3][0][0], &dLTrue[3][1][0], &dLTrue[3][2][0]},
254 {nullptr, nullptr, nullptr}};
255
256 for (unsigned int f = 0; f < 2; ++f)
257 {
258 for (unsigned int d = 0; d < 2; ++d)
259 for (int i = -1; i < mesh->n[f][d] + 1; ++i)
260 ASSERT_NEAR(dL[f][d][i], mesh->dL[f][d][i], 1e-12)
261 << "dL incorrect at field " << f << ", direction " << d
262 << ", index " << i << std::endl;
263
264 ASSERT_NEAR(dL[f][2][0], mesh->dL[f][2][0], 1e-12)
265 << "dL incorrect at field " << f << ", direction " << 2
266 << ", index " << 0 << std::endl;
267 }
268
269 for (unsigned int f = 2; f < 4; ++f)
270 for (unsigned int d = 0; d < 3; ++d)
271 for (int i = 0; i < mesh->n[f][d]; ++i)
272 ASSERT_NEAR(dL[f][d][i], mesh->dL[f][d][i], 1e-12)
273 << "dL incorrect at field " << f << ", direction " << d
274 << ", index " << i << std::endl;
275
276 for (unsigned int d = 0; d < 3; ++d)
277 ASSERT_FALSE(mesh->dL[4][d])
278 << "The dL of the vertex mesh is not nullptr." << std::endl;
279}
280
281// check UN
282TEST_F(CartesainMeshTest2D_YPeriodic, check_UN)
283{
284 ASSERT_EQ(253, mesh->UN)
285 << "The total number of velocity points is not correct." << std::endl;
286}
287
288// check pN
289TEST_F(CartesainMeshTest2D_YPeriodic, check_pN)
290{
291 ASSERT_EQ(132, mesh->pN)
292 << "The total number of pressure points is not correct." << std::endl;
293}
294
295// check da
296TEST_F(CartesainMeshTest2D_YPeriodic, check_da)
297{
298 PetscErrorCode ierr;
299 DMType type;
300
301 ierr = DMGetType(mesh->da[0], &type);
302 ASSERT_FALSE(ierr);
303 ASSERT_STREQ(type, DMDA)
304 << "u-velocity mesh is " << type << ", not expected da" << std::endl;
305
306 ierr = DMGetType(mesh->da[1], &type);
307 ASSERT_FALSE(ierr);
308 ASSERT_STREQ(type, DMDA)
309 << "v-velocity mesh is " << type << ", not expected da" << std::endl;
310
311 ierr = DMGetType(mesh->da[3], &type);
312 ASSERT_FALSE(ierr);
313 ASSERT_STREQ(type, DMDA)
314 << "v-velocity mesh is " << type << ", not expected da" << std::endl;
315
316 ASSERT_FALSE(mesh->da[2])
317 << "w-velocity mesh should not be initialized." << std::endl;
318 ASSERT_FALSE(mesh->da[4])
319 << "vertex mesh should not be initialized." << std::endl;
320}
321
322// check nProc
323TEST_F(CartesainMeshTest2D_YPeriodic, check_nProc)
324{
325 PetscErrorCode ierr;
326
327 PetscInt nProc[3];
328
329 ierr = DMDAGetInfo(da[3], nullptr, nullptr, nullptr, nullptr, &nProc[0],
330 &nProc[1], &nProc[2], nullptr, nullptr, nullptr, nullptr,
331 nullptr, nullptr);
332 ASSERT_FALSE(ierr);
333
334 ASSERT_EQ(nProc[0], mesh->nProc[0]);
335 ASSERT_EQ(nProc[1], mesh->nProc[1]);
336 ASSERT_EQ(nProc[2], mesh->nProc[2]);
337}
338
339// check bg
340TEST_F(CartesainMeshTest2D_YPeriodic, check_bg)
341{
342 PetscErrorCode ierr;
343 for (unsigned int f = 0; f < 5; f++)
344 {
345 if ((f == 2) || (f == 4)) continue;
346 DMDALocalInfo info;
347 ierr = DMDAGetLocalInfo(da[f], &info);
348 ASSERT_FALSE(ierr);
349 ASSERT_EQ(info.xs, mesh->bg[f][0]) << "f = " << f << std::endl;
350 ASSERT_EQ(info.ys, mesh->bg[f][1]) << "f = " << f << std::endl;
351 ASSERT_EQ(info.zs, mesh->bg[f][2]) << "f = " << f << std::endl;
352 }
353
354 ASSERT_EQ(0, mesh->bg[2][0]) << "f = " << 2 << std::endl;
355 ASSERT_EQ(0, mesh->bg[2][1]) << "f = " << 2 << std::endl;
356 ASSERT_EQ(0, mesh->bg[2][2]) << "f = " << 2 << std::endl;
357}
358
359// check ed
360TEST_F(CartesainMeshTest2D_YPeriodic, check_ed)
361{
362 PetscErrorCode ierr;
363 for (unsigned int f = 0; f < 5; f++)
364 {
365 if ((f == 2) || (f == 4)) continue;
366 DMDALocalInfo info;
367 ierr = DMDAGetLocalInfo(da[f], &info);
368 ASSERT_FALSE(ierr);
369 ASSERT_EQ(info.xs + info.xm, mesh->ed[f][0])
370 << "f = " << f << std::endl;
371 ASSERT_EQ(info.ys + info.ym, mesh->ed[f][1])
372 << "f = " << f << std::endl;
373 ASSERT_EQ(info.zs + info.zm, mesh->ed[f][2])
374 << "f = " << f << std::endl;
375 }
376
377 ASSERT_EQ(1, mesh->ed[2][0]) << "f = " << 2 << std::endl;
378 ASSERT_EQ(1, mesh->ed[2][1]) << "f = " << 2 << std::endl;
379 ASSERT_EQ(1, mesh->ed[2][2]) << "f = " << 2 << std::endl;
380}
381
382// check m
383TEST_F(CartesainMeshTest2D_YPeriodic, check_m)
384{
385 PetscErrorCode ierr;
386 for (unsigned int f = 0; f < 5; f++)
387 {
388 if ((f == 2) || (f == 4)) continue;
389 DMDALocalInfo info;
390 ierr = DMDAGetLocalInfo(da[f], &info);
391 ASSERT_FALSE(ierr);
392 ASSERT_EQ(info.xm, mesh->m[f][0]) << "f = " << f << std::endl;
393 ASSERT_EQ(info.ym, mesh->m[f][1]) << "f = " << f << std::endl;
394 ASSERT_EQ(info.zm, mesh->m[f][2]) << "f = " << f << std::endl;
395 }
396
397 ASSERT_EQ(0, mesh->m[2][0]) << "f = " << 2 << std::endl;
398 ASSERT_EQ(0, mesh->m[2][1]) << "f = " << 2 << std::endl;
399 ASSERT_EQ(0, mesh->m[2][2]) << "f = " << 2 << std::endl;
400}
401
402// check UNLocal
403TEST_F(CartesainMeshTest2D_YPeriodic, check_UNLocal)
404{
405 int N = mesh->m[0][0] * mesh->m[0][1] + mesh->m[1][0] * mesh->m[1][1];
406 ASSERT_EQ(N, mesh->UNLocal);
407}
408
409// check pNLocal
410TEST_F(CartesainMeshTest2D_YPeriodic, check_pNLocal)
411{
412 int N = mesh->m[3][0] * mesh->m[3][1];
413 ASSERT_EQ(N, mesh->pNLocal);
414}
415
416// check UPack
417TEST_F(CartesainMeshTest2D_YPeriodic, check_UPack)
418{
419 PetscErrorCode ierr;
420 DMType type;
421
422 ierr = DMGetType(mesh->UPack, &type);
423 ASSERT_FALSE(ierr);
424 ASSERT_STREQ(type, DMCOMPOSITE);
425
426 DM dms[3]{nullptr, nullptr, nullptr};
427 ierr = DMCompositeGetEntriesArray(mesh->UPack, dms);
428 ASSERT_FALSE(ierr);
429 ASSERT_EQ(dms[0], mesh->da[0]);
430 ASSERT_EQ(dms[1], mesh->da[1]);
431 ASSERT_EQ(dms[2], mesh->da[2]);
432
433 ierr = DMDestroy(&dms[0]);
434 ASSERT_FALSE(ierr);
435 ierr = DMDestroy(&dms[1]);
436 ASSERT_FALSE(ierr);
437}
438
439// check mpi
440TEST_F(CartesainMeshTest2D_YPeriodic, check_mpi)
441{
442 PetscErrorCode ierr;
443 PetscMPIInt rank, size;
444
445 ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);
446 ASSERT_FALSE(ierr);
447 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
448 ASSERT_FALSE(ierr);
449
450 ASSERT_EQ(PETSC_COMM_WORLD, mesh->comm);
451 ASSERT_EQ(size, mesh->mpiSize);
452 ASSERT_EQ(rank, mesh->mpiRank);
453}
454
455// check getLocalIndex with i, j, k
456TEST_F(CartesainMeshTest2D_YPeriodic, check_getLocalIndex_ijk)
457{
458 for (int f = 0; f < 5; ++f)
459 {
460 switch (f)
461 {
462 case 2:
463 case 4:
464 break;
465 default:
466 for (int k = mesh->bg[f][2]; k < mesh->ed[f][2]; k++)
467 for (int j = mesh->bg[f][1] - 1; j < mesh->ed[f][1] + 1;
468 j++)
469 for (int i = mesh->bg[f][0] - 1; i < mesh->ed[f][0] + 1;
470 i++)
471 {
472 PetscInt idx;
473 mesh->getLocalIndex(f, i, j, k, idx);
474
475 PetscInt i2, j2, k2, idx2;
476
477 i2 = i - mesh->bg[f][0] + 1;
478 j2 = j - mesh->bg[f][1] + 1;
479 k2 = k - mesh->bg[f][2];
480
481 idx2 = i2 + j2 * (mesh->m[f][0] + 2) +
482 k2 * ((mesh->m[f][0] + 2) *
483 (mesh->m[f][1] + 2));
484
485 ASSERT_EQ(idx2, idx)
486 << "getLocalIndex wrong at field " << f
487 << " and (" << i << ", " << j << ", " << k
488 << ") " << std::endl;
489 }
490 }
491 }
492}
493
494// check getNaturalIndex with i, j, k
495TEST_F(CartesainMeshTest2D_YPeriodic, check_getNaturalIndex_ijk)
496{
497 for (int f = 0; f < 5; ++f)
498 {
499 // internal points
500 for (int k = 0; k < mesh->n[f][2]; k++)
501 for (int j = 0; j < mesh->n[f][1]; j++)
502 for (int i = 0; i < mesh->n[f][0]; i++)
503 {
504 PetscInt idx, idx2;
505 mesh->getNaturalIndex(f, i, j, k, idx);
506
507 idx2 = i + j * mesh->n[f][0] +
508 k * (mesh->n[f][0] * mesh->n[f][1]);
509
510 ASSERT_EQ(idx2, idx)
511 << "getNatrualIndex wrong at field " << f
512 << " and internal points (" << i << ", " << j << ", "
513 << k << ") " << std::endl;
514 }
515
516 // ghost points at top and bottom
517 for (int k = 0; k < mesh->n[f][2]; k++)
518 for (int i = 0; i < mesh->n[f][0]; i++)
519 {
520 PetscInt idx, idx2, j;
521
522 j = mesh->n[f][1];
523 mesh->getNaturalIndex(f, i, j, k, idx);
524 idx2 = i + k * (mesh->n[f][0] * mesh->n[f][1]);
525 ASSERT_EQ(idx2, idx) << "getNatrualIndex wrong at field " << f
526 << " and top ghost point (" << i << ", "
527 << j << ", " << k << ") " << std::endl;
528
529 j = -1;
530 mesh->getNaturalIndex(f, i, j, k, idx);
531 idx2 = i + (mesh->n[f][1] - 1) * mesh->n[f][0] +
532 k * (mesh->n[f][0] * mesh->n[f][1]);
533 ASSERT_EQ(idx2, idx) << "getNatrualIndex wrong at field " << f
534 << " and bottom ghost point (" << i << ", "
535 << j << ", " << k << ") " << std::endl;
536 }
537
538 // ghost points at right and left
539 for (int k = 0; k < mesh->n[f][2]; k++)
540 for (int j = 0; j < mesh->n[f][1]; j++)
541 {
542 PetscInt idx, i;
543
544 i = mesh->n[f][0];
545 mesh->getNaturalIndex(f, i, j, k, idx);
546 ASSERT_EQ(-1, idx) << "getNatrualIndex wrong at field " << f
547 << " and right ghost point (" << i << ", "
548 << j << ", " << k << ") " << std::endl;
549
550 i = -1;
551 mesh->getNaturalIndex(f, i, j, k, idx);
552 ASSERT_EQ(-1, idx) << "getNatrualIndex wrong at field " << f
553 << " and left ghost point (" << i << ", "
554 << j << ", " << k << ") " << std::endl;
555 }
556
557 // corners
558 for (int k = 0; k < mesh->n[f][2]; k++)
559 {
560 PetscInt idx;
561
562 mesh->getNaturalIndex(f, -1, -1, 0, idx);
563 ASSERT_EQ(-1, idx) << "getNatrualIndex wrong at field " << f
564 << " and corner ghost point (" << -1 << ", "
565 << -1 << ", " << k << ") " << std::endl;
566
567 mesh->getNaturalIndex(f, -1, mesh->n[f][1], 0, idx);
568 ASSERT_EQ(-1, idx)
569 << "getNatrualIndex wrong at field " << f
570 << " and corner ghost point (" << -1 << ", " << mesh->n[f][1]
571 << ", " << k << ") " << std::endl;
572
573 mesh->getNaturalIndex(f, mesh->n[f][0], -1, 0, idx);
574 ASSERT_EQ(-1, idx) << "getNatrualIndex wrong at field " << f
575 << " and corner ghost point (" << mesh->n[f][0]
576 << ", " << -1 << ", " << k << ") " << std::endl;
577
578 mesh->getNaturalIndex(f, mesh->n[f][0], mesh->n[f][1], 0, idx);
579 ASSERT_EQ(-1, idx)
580 << "getNatrualIndex wrong at field " << f
581 << " and corner ghost point (" << mesh->n[f][0] << ", "
582 << mesh->n[f][1] << ", " << k << ") " << std::endl;
583 }
584 }
585}
586
587// check getGlobalIndex with i, j, k of internal points
588TEST_F(CartesainMeshTest2D_YPeriodic, check_getGlobalIndex_ijk_internal)
589{
590 for (int f = 0; f < 5; ++f)
591 {
592 if ((f == 2) || (f == 4)) continue;
593
594 PetscErrorCode ierr;
595
596 ISLocalToGlobalMapping mapping;
597 ierr = DMGetLocalToGlobalMapping(mesh->da[f], &mapping);
598 ASSERT_FALSE(ierr);
599
600 Vec vini, vinj, vink;
601 ierr = DMGetGlobalVector(mesh->da[f], &vini);
602 ASSERT_FALSE(ierr);
603 ierr = DMGetGlobalVector(mesh->da[f], &vinj);
604 ASSERT_FALSE(ierr);
605 ierr = DMGetGlobalVector(mesh->da[f], &vink);
606 ASSERT_FALSE(ierr);
607
608 for (int k = mesh->bg[f][2]; k < mesh->ed[f][2]; k++)
609 for (int j = mesh->bg[f][1]; j < mesh->ed[f][1]; j++)
610 for (int i = mesh->bg[f][0]; i < mesh->ed[f][0]; i++)
611 {
612 PetscInt idx;
613
614 ierr = DMDAConvertToCell(mesh->da[f], {k, j, i, 0}, &idx);
615 ASSERT_FALSE(ierr);
616 ierr = ISLocalToGlobalMappingApply(mapping, 1, &idx, &idx);
617 ASSERT_FALSE(ierr);
618
619 ierr = VecSetValue(vini, idx, i, INSERT_VALUES);
620 ASSERT_FALSE(ierr);
621 ierr = VecSetValue(vinj, idx, j, INSERT_VALUES);
622 ASSERT_FALSE(ierr);
623 ierr = VecSetValue(vink, idx, k, INSERT_VALUES);
624 ASSERT_FALSE(ierr);
625 }
626 ierr = VecAssemblyBegin(vini);
627 ASSERT_FALSE(ierr);
628 ierr = VecAssemblyBegin(vinj);
629 ASSERT_FALSE(ierr);
630 ierr = VecAssemblyBegin(vink);
631 ASSERT_FALSE(ierr);
632 ierr = VecAssemblyEnd(vini);
633 ASSERT_FALSE(ierr);
634 ierr = VecAssemblyEnd(vinj);
635 ASSERT_FALSE(ierr);
636 ierr = VecAssemblyEnd(vink);
637 ASSERT_FALSE(ierr);
638
639 Vec vouti, voutj, voutk;
640 VecScatter ctx;
641 ierr = VecScatterCreateToAll(vini, &ctx, &vouti);
642 ASSERT_FALSE(ierr);
643 ierr = VecDuplicate(vouti, &voutj);
644 ASSERT_FALSE(ierr);
645 ierr = VecDuplicate(vouti, &voutk);
646 ASSERT_FALSE(ierr);
647 ierr =
648 VecScatterBegin(ctx, vini, vouti, INSERT_VALUES, SCATTER_FORWARD);
649 ASSERT_FALSE(ierr);
650 ierr = VecScatterEnd(ctx, vini, vouti, INSERT_VALUES, SCATTER_FORWARD);
651 ASSERT_FALSE(ierr);
652 ierr =
653 VecScatterBegin(ctx, vinj, voutj, INSERT_VALUES, SCATTER_FORWARD);
654 ASSERT_FALSE(ierr);
655 ierr = VecScatterEnd(ctx, vinj, voutj, INSERT_VALUES, SCATTER_FORWARD);
656 ASSERT_FALSE(ierr);
657 ierr =
658 VecScatterBegin(ctx, vink, voutk, INSERT_VALUES, SCATTER_FORWARD);
659 ASSERT_FALSE(ierr);
660 ierr = VecScatterEnd(ctx, vink, voutk, INSERT_VALUES, SCATTER_FORWARD);
661 ASSERT_FALSE(ierr);
662
663 // internal points
664 for (int k = 0; k < mesh->n[f][2]; k++)
665 for (int j = 0; j < mesh->n[f][1]; j++)
666 for (int i = 0; i < mesh->n[f][0]; i++)
667 {
668 PetscInt idx;
669 ierr = mesh->getGlobalIndex(f, i, j, k, idx);
670 ASSERT_FALSE(ierr);
671
672 PetscReal vi, vj, vk;
673 ierr = VecGetValues(vouti, 1, &idx, &vi);
674 ASSERT_FALSE(ierr);
675 ierr = VecGetValues(voutj, 1, &idx, &vj);
676 ASSERT_FALSE(ierr);
677 ierr = VecGetValues(voutk, 1, &idx, &vk);
678 ASSERT_FALSE(ierr);
679
680 ASSERT_EQ(vi, PetscReal(i));
681 ASSERT_EQ(vj, PetscReal(j));
682 ASSERT_EQ(vk, PetscReal(k));
683 }
684
685 ierr = VecScatterDestroy(&ctx);
686 ASSERT_FALSE(ierr);
687 ierr = VecDestroy(&voutk);
688 ASSERT_FALSE(ierr);
689 ierr = VecDestroy(&voutj);
690 ASSERT_FALSE(ierr);
691 ierr = VecDestroy(&vouti);
692 ASSERT_FALSE(ierr);
693 ierr = VecDestroy(&vink);
694 ASSERT_FALSE(ierr);
695 ierr = VecDestroy(&vinj);
696 ASSERT_FALSE(ierr);
697 ierr = VecDestroy(&vini);
698 ASSERT_FALSE(ierr);
699
700 // ierr = ISLocalToGlobalMappingDestroy(&mapping); ASSERT_FALSE(ierr);
701 // Seem there's a bug in ISLocalToGlobalMappingDestroy. It should just
702 // decrease the reference number by one. Instead, it now destroys the
703 // underlying object completely....
704 mapping = PETSC_NULL;
705 }
706}
707
708// check getGlobalIndex with i, j, k of ghost points
709TEST_F(CartesainMeshTest2D_YPeriodic, check_getGlobalIndex_ijk_ghosts)
710{
711 for (int f = 0; f < 5; ++f)
712 {
713 if ((f == 2) || (f == 4)) continue;
714
715 // ghost points at top and bottom
716 for (int k = 0; k < mesh->n[f][2]; k++)
717 for (int i = 0; i < mesh->n[f][0]; i++)
718 {
719 PetscInt idx, idx2, j;
720
721 j = mesh->n[f][1];
722 mesh->getGlobalIndex(f, i, j, k, idx);
723
724 // internal part should have already passed
725 mesh->getGlobalIndex(f, i, 0, k, idx2);
726 ASSERT_EQ(idx2, idx)
727 << "getGlobalIndex wrong at field " << f << " and (" << i
728 << ", " << j << ", " << k << ") " << std::endl;
729
730 j = -1;
731 mesh->getGlobalIndex(f, i, j, k, idx);
732
733 // internal part should have already passed
734 mesh->getGlobalIndex(f, i, mesh->n[f][1] - 1, k, idx2);
735 ASSERT_EQ(idx2, idx)
736 << "getGlobalIndex wrong at field " << f << " and (" << i
737 << ", " << j << ", " << k << ") " << std::endl;
738 }
739
740 // ghost points at right and left
741 for (int k = 0; k < mesh->n[f][2]; k++)
742 for (int j = 0; j < mesh->n[f][1]; j++)
743 {
744 PetscInt idx, i;
745
746 i = mesh->n[f][0];
747 mesh->getGlobalIndex(f, i, j, k, idx);
748 ASSERT_EQ(-1, idx)
749 << "getGlobalIndex wrong at field " << f << " and (" << i
750 << ", " << j << ", " << k << ") " << std::endl;
751
752 i = -1;
753 mesh->getGlobalIndex(f, i, j, k, idx);
754 ASSERT_EQ(-1, idx)
755 << "getGlobalIndex wrong at field " << f << " and (" << i
756 << ", " << j << ", " << k << ") " << std::endl;
757 }
758
759 // corners
760 for (int k = 0; k < mesh->n[f][2]; k++)
761 {
762 PetscInt idx;
763
764 mesh->getGlobalIndex(f, -1, -1, 0, idx);
765 ASSERT_EQ(-1, idx)
766 << "getGlobalIndex wrong at field " << f << " and (" << -1
767 << ", " << -1 << ", " << k << ") " << std::endl;
768
769 mesh->getGlobalIndex(f, -1, mesh->n[f][1], 0, idx);
770 ASSERT_EQ(-1, idx)
771 << "getGlobalIndex wrong at field " << f << " and (" << -1
772 << ", " << mesh->n[f][1] << ", " << k << ") " << std::endl;
773
774 mesh->getGlobalIndex(f, mesh->n[f][0], -1, 0, idx);
775 ASSERT_EQ(-1, idx) << "getGlobalIndex wrong at field " << f
776 << " and (" << mesh->n[f][0] << ", " << -1
777 << ", " << k << ") " << std::endl;
778
779 mesh->getGlobalIndex(f, mesh->n[f][0], mesh->n[f][1], 0, idx);
780 ASSERT_EQ(-1, idx)
781 << "getGlobalIndex wrong at field " << f << " and ("
782 << mesh->n[f][0] << ", " << mesh->n[f][1] << ", " << k << ") "
783 << std::endl;
784 }
785 }
786}
787
788// check getPackedGlobalIndex with i, j, k of internal points
789TEST_F(CartesainMeshTest2D_YPeriodic, check_getPackedGlobalIndex_ijk_internal)
790{
791 PetscErrorCode ierr;
792
793 DM UPack;
794
795 ierr = DMCompositeCreate(PETSC_COMM_WORLD, &UPack);
796 ASSERT_FALSE(ierr);
797 ierr = DMCompositeAddDM(UPack, da[0]);
798 ASSERT_FALSE(ierr);
799 ierr = DMCompositeAddDM(UPack, da[1]);
800 ASSERT_FALSE(ierr);
801 ierr = DMSetUp(UPack);
802 ASSERT_FALSE(ierr);
803
804 Vec fp, ip, jp, kp;
805 ierr = DMCreateGlobalVector(mesh->UPack, &fp);
806 ASSERT_FALSE(ierr);
807 ierr = DMCreateGlobalVector(mesh->UPack, &ip);
808 ASSERT_FALSE(ierr);
809 ierr = DMCreateGlobalVector(mesh->UPack, &jp);
810 ASSERT_FALSE(ierr);
811 ierr = DMCreateGlobalVector(mesh->UPack, &kp);
812 ASSERT_FALSE(ierr);
813
814 Vec fv[2], iv[2], jv[2], kv[2];
815 ierr = DMCompositeGetAccessArray(mesh->UPack, fp, 2, nullptr, fv);
816 ASSERT_FALSE(ierr);
817 ierr = DMCompositeGetAccessArray(mesh->UPack, ip, 2, nullptr, iv);
818 ASSERT_FALSE(ierr);
819 ierr = DMCompositeGetAccessArray(mesh->UPack, jp, 2, nullptr, jv);
820 ASSERT_FALSE(ierr);
821 ierr = DMCompositeGetAccessArray(mesh->UPack, kp, 2, nullptr, kv);
822 ASSERT_FALSE(ierr);
823
824 for (int f = 0; f < 2; ++f)
825 {
826 ISLocalToGlobalMapping mapping;
827 ierr = DMGetLocalToGlobalMapping(mesh->da[f], &mapping);
828 ASSERT_FALSE(ierr);
829
830 for (int k = mesh->bg[f][2]; k < mesh->ed[f][2]; ++k)
831 for (int j = mesh->bg[f][1]; j < mesh->ed[f][1]; ++j)
832 for (int i = mesh->bg[f][0]; i < mesh->ed[f][0]; ++i)
833 {
834 PetscInt idx;
835 ierr = DMDAConvertToCell(mesh->da[f], {k, j, i, 0}, &idx);
836 ASSERT_FALSE(ierr);
837 ierr = ISLocalToGlobalMappingApply(mapping, 1, &idx, &idx);
838 ASSERT_FALSE(ierr);
839
840 ierr = VecSetValue(fv[f], idx, (PetscReal)f, INSERT_VALUES);
841 ASSERT_FALSE(ierr);
842 ierr = VecSetValue(iv[f], idx, (PetscReal)i, INSERT_VALUES);
843 ASSERT_FALSE(ierr);
844 ierr = VecSetValue(jv[f], idx, (PetscReal)j, INSERT_VALUES);
845 ASSERT_FALSE(ierr);
846 ierr = VecSetValue(kv[f], idx, (PetscReal)k, INSERT_VALUES);
847 ASSERT_FALSE(ierr);
848 }
849
850 ierr = VecAssemblyBegin(fv[f]);
851 ASSERT_FALSE(ierr);
852 ierr = VecAssemblyBegin(iv[f]);
853 ASSERT_FALSE(ierr);
854 ierr = VecAssemblyBegin(jv[f]);
855 ASSERT_FALSE(ierr);
856 ierr = VecAssemblyBegin(kv[f]);
857 ASSERT_FALSE(ierr);
858 ierr = VecAssemblyEnd(fv[f]);
859 ASSERT_FALSE(ierr);
860 ierr = VecAssemblyEnd(iv[f]);
861 ASSERT_FALSE(ierr);
862 ierr = VecAssemblyEnd(jv[f]);
863 ASSERT_FALSE(ierr);
864 ierr = VecAssemblyEnd(kv[f]);
865 ASSERT_FALSE(ierr);
866
867 ierr = ISLocalToGlobalMappingDestroy(&mapping);
868 ASSERT_FALSE(ierr);
869 }
870
871 ierr = DMCompositeRestoreAccessArray(mesh->UPack, fp, 2, nullptr, fv);
872 ASSERT_FALSE(ierr);
873 ierr = DMCompositeRestoreAccessArray(mesh->UPack, ip, 2, nullptr, iv);
874 ASSERT_FALSE(ierr);
875 ierr = DMCompositeRestoreAccessArray(mesh->UPack, jp, 2, nullptr, jv);
876 ASSERT_FALSE(ierr);
877 ierr = DMCompositeRestoreAccessArray(mesh->UPack, kp, 2, nullptr, kv);
878 ASSERT_FALSE(ierr);
879
880 VecScatter ctx;
881 Vec fps, ips, jps, kps;
882
883 ierr = VecScatterCreateToAll(fp, &ctx, &fps);
884 ASSERT_FALSE(ierr);
885 ierr = VecDuplicate(fps, &ips);
886 ASSERT_FALSE(ierr);
887 ierr = VecDuplicate(fps, &jps);
888 ASSERT_FALSE(ierr);
889 ierr = VecDuplicate(fps, &kps);
890 ASSERT_FALSE(ierr);
891
892 ierr = VecScatterBegin(ctx, fp, fps, INSERT_VALUES, SCATTER_FORWARD);
893 ASSERT_FALSE(ierr);
894 ierr = VecScatterEnd(ctx, fp, fps, INSERT_VALUES, SCATTER_FORWARD);
895 ASSERT_FALSE(ierr);
896 ierr = VecScatterBegin(ctx, ip, ips, INSERT_VALUES, SCATTER_FORWARD);
897 ASSERT_FALSE(ierr);
898 ierr = VecScatterEnd(ctx, ip, ips, INSERT_VALUES, SCATTER_FORWARD);
899 ASSERT_FALSE(ierr);
900 ierr = VecScatterBegin(ctx, jp, jps, INSERT_VALUES, SCATTER_FORWARD);
901 ASSERT_FALSE(ierr);
902 ierr = VecScatterEnd(ctx, jp, jps, INSERT_VALUES, SCATTER_FORWARD);
903 ASSERT_FALSE(ierr);
904 ierr = VecScatterBegin(ctx, kp, kps, INSERT_VALUES, SCATTER_FORWARD);
905 ASSERT_FALSE(ierr);
906 ierr = VecScatterEnd(ctx, kp, kps, INSERT_VALUES, SCATTER_FORWARD);
907 ASSERT_FALSE(ierr);
908
909 ierr = VecScatterDestroy(&ctx);
910 ASSERT_FALSE(ierr);
911 ierr = VecDestroy(&kp);
912 ASSERT_FALSE(ierr);
913 ierr = VecDestroy(&jp);
914 ASSERT_FALSE(ierr);
915 ierr = VecDestroy(&ip);
916 ASSERT_FALSE(ierr);
917 ierr = VecDestroy(&fp);
918 ASSERT_FALSE(ierr);
919
920 for (int f = 0; f < 2; ++f)
921 {
922 // internal points
923 for (int k = 0; k < mesh->n[f][2]; k++)
924 for (int j = 0; j < mesh->n[f][1]; j++)
925 for (int i = 0; i < mesh->n[f][0]; i++)
926 {
927 PetscInt idx;
928 ierr = mesh->getPackedGlobalIndex(f, i, j, k, idx);
929 ASSERT_FALSE(ierr);
930
931 PetscReal vf, vi, vj, vk;
932 ierr = VecGetValues(fps, 1, &idx, &vf);
933 ASSERT_FALSE(ierr);
934 ierr = VecGetValues(ips, 1, &idx, &vi);
935 ASSERT_FALSE(ierr);
936 ierr = VecGetValues(jps, 1, &idx, &vj);
937 ASSERT_FALSE(ierr);
938 ierr = VecGetValues(kps, 1, &idx, &vk);
939 ASSERT_FALSE(ierr);
940
941 ASSERT_EQ(vf, f) << f << ".(" << i << ", " << j << ", " << k
942 << ")" << std::endl;
943 ASSERT_EQ(vi, i) << f << ".(" << i << ", " << j << ", " << k
944 << ")" << std::endl;
945 ASSERT_EQ(vj, j) << f << ".(" << i << ", " << j << ", " << k
946 << ")" << std::endl;
947 ASSERT_EQ(vk, k) << f << ".(" << i << ", " << j << ", " << k
948 << ")" << std::endl;
949 }
950 }
951
952 ierr = VecDestroy(&kps);
953 ASSERT_FALSE(ierr);
954 ierr = VecDestroy(&jps);
955 ASSERT_FALSE(ierr);
956 ierr = VecDestroy(&ips);
957 ASSERT_FALSE(ierr);
958 ierr = VecDestroy(&fps);
959 ASSERT_FALSE(ierr);
960}
961
962// check getPackedGlobalIndex with i, j, k of ghosts
963TEST_F(CartesainMeshTest2D_YPeriodic, check_getPackedGlobalIndex_ijk_ghost)
964{
965 for (int f = 0; f < 2; ++f)
966 {
967 // ghost points at top and bottom
968 for (int k = 0; k < mesh->n[f][2]; k++)
969 for (int i = 0; i < mesh->n[f][0]; i++)
970 {
971 PetscInt idx, idx2, j;
972
973 j = mesh->n[f][1];
974 mesh->getPackedGlobalIndex(f, i, j, k, idx);
975
976 // internal part should have already passed
977 mesh->getPackedGlobalIndex(f, i, 0, k, idx2);
978 ASSERT_EQ(idx2, idx)
979 << "getPackedGlobalIndex wrong at field " << f << " and ("
980 << i << ", " << j << ", " << k << ") " << std::endl;
981
982 j = -1;
983 mesh->getPackedGlobalIndex(f, i, j, k, idx);
984
985 // internal part should have already passed
986 mesh->getPackedGlobalIndex(f, i, mesh->n[f][1] - 1, k, idx2);
987 ASSERT_EQ(idx2, idx)
988 << "getPackedGlobalIndex wrong at field " << f << " and ("
989 << i << ", " << j << ", " << k << ") " << std::endl;
990 }
991
992 // ghost points at right and left
993 for (int k = 0; k < mesh->n[f][2]; k++)
994 for (int j = 0; j < mesh->n[f][1]; j++)
995 {
996 PetscInt idx, i;
997
998 i = mesh->n[f][0];
999 mesh->getPackedGlobalIndex(f, i, j, k, idx);
1000 ASSERT_EQ(-1, idx)
1001 << "getPackedGlobalIndex wrong at field " << f << " and ("
1002 << i << ", " << j << ", " << k << ") " << std::endl;
1003
1004 i = -1;
1005 mesh->getPackedGlobalIndex(f, i, j, k, idx);
1006 ASSERT_EQ(-1, idx)
1007 << "getPackedGlobalIndex wrong at field " << f << " and ("
1008 << i << ", " << j << ", " << k << ") " << std::endl;
1009 }
1010
1011 // corners
1012 for (int k = 0; k < mesh->n[f][2]; k++)
1013 {
1014 PetscInt idx;
1015
1016 mesh->getPackedGlobalIndex(f, -1, -1, 0, idx);
1017 ASSERT_EQ(-1, idx)
1018 << "getPackedGlobalIndex wrong at field " << f << " and (" << -1
1019 << ", " << -1 << ", " << k << ") " << std::endl;
1020
1021 mesh->getPackedGlobalIndex(f, -1, mesh->n[f][1], 0, idx);
1022 ASSERT_EQ(-1, idx)
1023 << "getPackedGlobalIndex wrong at field " << f << " and (" << -1
1024 << ", " << mesh->n[f][1] << ", " << k << ") " << std::endl;
1025
1026 mesh->getPackedGlobalIndex(f, mesh->n[f][0], -1, 0, idx);
1027 ASSERT_EQ(-1, idx) << "getPackedGlobalIndex wrong at field " << f
1028 << " and (" << mesh->n[f][0] << ", " << -1
1029 << ", " << k << ") " << std::endl;
1030
1031 mesh->getPackedGlobalIndex(f, mesh->n[f][0], mesh->n[f][1], 0, idx);
1032 ASSERT_EQ(-1, idx)
1033 << "getPackedGlobalIndex wrong at field " << f << " and ("
1034 << mesh->n[f][0] << ", " << mesh->n[f][1] << ", " << k << ") "
1035 << std::endl;
1036 }
1037 }
1038}
TEST_F(CartesainMeshTest2D_YPeriodic, check_dim)
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< BoolVec1D > BoolVec2D
2D std::vector holding PetscBool.
Definition: type.h:149
std::vector< RealVec1D > RealVec2D
2D std::vector holding PetscReal.
Definition: type.h:142
std::vector< RealVec2D > RealVec3D
3D std::vector holding PetscReal.
Definition: type.h:144
Prototype of mesh::MeshBase, type::Mesh, and factory function.
Definition: parser.cpp:45
std::vector< GhostedVec2D > GhostedVec3D
a vector of vector pointers to mimic ghosted 2D vectors. type
Definition: type.h:157