Nektar++
OutputVtk.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: OutputVtk.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// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: VTK file format output using VTK library.
32//
33////////////////////////////////////////////////////////////////////////////////
34
35#include <boost/format.hpp>
36
37#include "OutputVtk.h"
40
41#include <vtkCellType.h>
42#include <vtkDoubleArray.h>
43#include <vtkInformation.h>
44#include <vtkMultiBlockDataSet.h>
45#include <vtkPointData.h>
46#include <vtkPoints.h>
47#include <vtkXMLMultiBlockDataWriter.h>
48#include <vtkXMLUnstructuredGridWriter.h>
49
50namespace Nektar::FieldUtils
51{
52
54 ModuleKey(eOutputModule, "vtu"), OutputVtk::create, "Writes a VTU file.");
55
57{
59 m_config["legacy"] =
60 ConfigOption(true, "0", "Output using legacy manual file writers");
61 m_config["highorder"] = ConfigOption(
62 true, "0", "Output using new high-order Lagrange elements");
63 m_config["multiblock"] = ConfigOption(
64 true, "0", "Output using multi-blocks to separate composites");
65 m_config["uncompress"] = ConfigOption(true, "0", "Uncompress xml sections");
66}
67
68// Anonymous namespace for spectral order -> VTK order mapping functions
69// be aware that these orderings are very similar to gmsh but not the same!
70namespace
71{
72// Hashing function for the ordering map taking two integers
73inline size_t key2(int i, int j)
74{
75 return (size_t)i << 32 | (unsigned int)j;
76}
77
78// Hashing function for the ordering map taking three integers
79inline size_t key3(int i, int j, int k)
80{
81 return (size_t)i << 10 ^ j << 5 ^ k;
82}
83
84// Map that takes (a,b,c) -> m
85typedef std::tuple<int, int, int> Mode;
86
87struct cmpop
88{
89 bool operator()(Mode const &a, Mode const &b) const
90 {
91 if (std::get<0>(a) < std::get<0>(b))
92 {
93 return true;
94 }
95 if (std::get<0>(a) > std::get<0>(b))
96 {
97 return false;
98 }
99 if (std::get<1>(a) < std::get<1>(b))
100 {
101 return true;
102 }
103 if (std::get<1>(a) > std::get<1>(b))
104 {
105 return false;
106 }
107 if (std::get<2>(a) < std::get<2>(b))
108 {
109 return true;
110 }
111
112 return false;
113 }
114};
115
116void Rotate(int nrot, std::vector<long long> &surfVerts)
117{
118 int n, i, j, cnt;
119 int np = static_cast<int>(
120 (sqrt(8.0 * static_cast<int>(surfVerts.size()) + 1.0) - 1) / 2);
121 std::vector<long long> tmp(np * np);
122
123 for (n = 0; n < nrot; ++n)
124 {
125 for (cnt = i = 0; i < np; ++i)
126 {
127 for (j = 0; j < np - i; ++j, cnt++)
128 {
129 tmp[i * np + j] = surfVerts[cnt];
130 }
131 }
132 for (cnt = i = 0; i < np; ++i)
133 {
134 for (j = 0; j < np - i; ++j, cnt++)
135 {
136 surfVerts[cnt] = tmp[(np - 1 - i - j) * np + i];
137 }
138 }
139 }
140}
141
142void Reflect(std::vector<long long> &surfVerts)
143{
144 int i, j, cnt;
145 int np = static_cast<int>(
146 (sqrt(8.0 * static_cast<double>(surfVerts.size()) + 1.0) - 1) / 2);
147 std::vector<long long> tmp(np * np);
148
149 for (cnt = i = 0; i < np; ++i)
150 {
151 for (j = 0; j < np - i; ++j, cnt++)
152 {
153 tmp[i * np + np - i - 1 - j] = surfVerts[cnt];
154 }
155 }
156
157 for (cnt = i = 0; i < np; ++i)
158 {
159 for (j = 0; j < np - i; ++j, cnt++)
160 {
161 surfVerts[cnt] = tmp[i * np + j];
162 }
163 }
164}
165
166void Align(std::vector<long long> thisVertId, std::vector<long long> vertId,
167 std::vector<long long> &surfVerts)
168{
169 if (vertId[0] == thisVertId[0])
170 {
171 if (vertId[1] == thisVertId[1] || vertId[1] == thisVertId[2])
172 {
173 if (vertId[1] == thisVertId[2])
174 {
175 Rotate(1, surfVerts);
176 Reflect(surfVerts);
177 }
178 }
179 }
180 else if (vertId[0] == thisVertId[1])
181 {
182 if (vertId[1] == thisVertId[0] || vertId[1] == thisVertId[2])
183 {
184 if (vertId[1] == thisVertId[0])
185 {
186 Reflect(surfVerts);
187 }
188 else
189 {
190 Rotate(2, surfVerts);
191 }
192 }
193 }
194 else if (vertId[0] == thisVertId[2])
195 {
196 if (vertId[1] == thisVertId[0] || vertId[1] == thisVertId[1])
197 {
198 if (vertId[1] == thisVertId[1])
199 {
200 Rotate(2, surfVerts);
201 Reflect(surfVerts);
202 }
203 else
204 {
205 Rotate(1, surfVerts);
206 }
207 }
208 }
209}
210
211std::vector<long long> quadTensorNodeOrdering(
212 const std::vector<long long> &nodes)
213{
214 int nN = static_cast<int>(nodes.size());
215 int n = static_cast<int>(sqrt(nN));
216
217 std::vector<long long> nodeList(nN);
218
219 // Vertices
220 nodeList[0] = nodes[0];
221 if (n > 1)
222 {
223 nodeList[n - 1] = nodes[1];
224 nodeList[n * n - 1] = nodes[2];
225 nodeList[n * (n - 1)] = nodes[3];
226 }
227
228 if (n > 2)
229 {
230 // Edge 0 -> 1
231 for (int i = 1; i < n - 1; ++i)
232 {
233 nodeList[i] = nodes[4 + i - 1];
234 }
235
236 // Edge 3 -> 2
237 for (int i = 1; i < n - 1; ++i)
238 {
239 nodeList[n * (n - 1) + i] = nodes[4 + 2 * (n - 2) + i - 1];
240 }
241
242 for (int j = 1; j < n - 1; ++j)
243 {
244 // Edge 0 -> 3
245 nodeList[n * (n - j - 1)] = nodes[4 + 3 * (n - 2) + n - 2 - j];
246
247 // Interior
248 for (int i = 1; i < n - 1; ++i)
249 {
250 nodeList[j * n + i] = nodes[(i - 3) - 2 * j + (3 + j) * n];
251 }
252
253 // Edge 1 -> 2
254 nodeList[(j + 1) * n - 1] = nodes[4 + (n - 2) + j - 1];
255 }
256 }
257
258 return nodeList;
259}
260
261std::vector<long long> triTensorNodeOrdering(
262 const std::vector<long long> &nodes)
263{
264 int nN = static_cast<int>(nodes.size());
265 int n = static_cast<int>(std::sqrt(2 * nN));
266
267 std::vector<long long> nodeList(nN);
268 int cnt2;
269
270 // Vertices
271 nodeList[0] = nodes[0];
272 if (n > 1)
273 {
274 nodeList[n - 1] = nodes[1];
275 nodeList[n * (n + 1) / 2 - 1] = nodes[2];
276 }
277
278 // Edges
279 int cnt = n;
280 for (int i = 1; i < n - 1; ++i)
281 {
282 nodeList[i] = nodes[3 + i - 1];
283 nodeList[cnt] = nodes[3 + 3 * (n - 2) - i];
284 nodeList[cnt + n - i - 1] = nodes[3 + (n - 2) + i - 1];
285 cnt += n - i;
286 }
287
288 // Interior (recursion)
289 if (n > 3)
290 {
291 // Reorder interior nodes
292 std::vector<long long> interior((n - 3) * (n - 2) / 2);
293 std::copy(nodes.begin() + 3 + 3 * (n - 2), nodes.end(),
294 interior.begin());
295 interior = triTensorNodeOrdering(interior);
296
297 // Copy into full node list
298 cnt = n;
299 cnt2 = 0;
300 for (int j = 1; j < n - 2; ++j)
301 {
302 for (int i = 0; i < n - j - 2; ++i)
303 {
304 nodeList[cnt + i + 1] = interior[cnt2 + i];
305 }
306 cnt += n - j;
307 cnt2 += n - 2 - j;
308 }
309 }
310
311 return nodeList;
312}
313
314std::vector<long long> tetTensorNodeOrdering(
315 const std::vector<long long> &nodes)
316{
317 int nN = static_cast<int>(nodes.size());
318 int n = static_cast<int>(
319 std::cbrt(3 * nN + std::sqrt(9 * nN * nN - 1 / 27)) +
320 std::cbrt(3 * nN - std::sqrt(9 * nN * nN - 1 / 27)) - 0.5);
321
322 std::vector<long long> nodeList(nN);
323 int nTri = n * (n + 1) / 2;
324 int nTet = n * (n + 1) * (n + 2) / 6;
325
326 // Vertices
327 nodeList[0] = nodes[0];
328 if (n == 1)
329 {
330 return nodeList;
331 }
332
333 nodeList[n - 1] = nodes[1];
334 nodeList[nTri - 1] = nodes[2];
335 nodeList[nTet - 1] = nodes[3];
336
337 if (n == 2)
338 {
339 return nodeList;
340 }
341
342 // Set up a map that takes (a,b,c) -> m to help us figure out where things
343 // are inside the tetrahedron.
344 std::map<Mode, int, cmpop> tmp;
345 for (int k = 0, cnt = 0; k < n; ++k)
346 {
347 for (int j = 0; j < n - k; ++j)
348 {
349 for (int i = 0; i < n - k - j; ++i)
350 {
351 tmp[Mode(i, j, k)] = cnt++;
352 }
353 }
354 }
355
356 // Edges first 3
357 for (int i = 1; i < n - 1; ++i)
358 {
359 int eI = i - 1;
360 nodeList[tmp[Mode(i, 0, 0)]] = nodes[4 + eI];
361 nodeList[tmp[Mode(n - 1 - i, i, 0)]] = nodes[4 + (n - 2) + eI];
362 nodeList[tmp[Mode(0, n - 1 - i, 0)]] = nodes[4 + 2 * (n - 2) + eI];
363 }
364
365 // Edges last 3 reversed (compared with NekMesh)
366 for (int i = 1; i < n - 1; ++i)
367 {
368 int eI = (n - 1 - i) - 1;
369 nodeList[tmp[Mode(0, 0, n - 1 - i)]] = nodes[4 + 3 * (n - 2) + eI];
370 nodeList[tmp[Mode(i, 0, n - 1 - i)]] = nodes[4 + 4 * (n - 2) + eI];
371 nodeList[tmp[Mode(0, i, n - 1 - i)]] = nodes[4 + 5 * (n - 2) + eI];
372 }
373
374 if (n == 3)
375 {
376 return nodeList;
377 }
378
379 // For faces, we use the triTensorNodeOrdering routine to make our lives
380 // slightly easier.
381 int nFacePts = (n - 3) * (n - 2) / 2;
382
383 // Grab face points and reorder into a tensor-product type format
384 std::vector<std::vector<long long>> tmpNodes(4);
385 int offset = 4 + 6 * (n - 2);
386
387 for (int i = 0; i < 4; ++i)
388 {
389 tmpNodes[i].resize(nFacePts);
390 for (int j = 0; j < nFacePts; ++j)
391 {
392 tmpNodes[i][j] = nodes[offset++];
393 }
394 tmpNodes[i] = triTensorNodeOrdering(tmpNodes[i]);
395 }
396
397 if (n > 4)
398 {
399 // Now align faces (different to NekMesh)
400 std::vector<long long> triVertId(3), toAlign(3);
401 triVertId[0] = 0;
402 triVertId[1] = 1;
403 triVertId[2] = 2;
404
405 toAlign[0] = 0;
406 toAlign[1] = 2;
407 toAlign[2] = 1;
408 Align(triVertId, toAlign, tmpNodes[2]);
409 Align(triVertId, toAlign, tmpNodes[3]);
410
411 toAlign[0] = 2;
412 toAlign[1] = 0;
413 toAlign[2] = 1;
414 Align(triVertId, toAlign, tmpNodes[1]);
415 }
416
417 // Reordered from NekMesh to put base last
418 for (int j = 1, cnt = 0; j < n - 2; ++j)
419 {
420 for (int i = 1; i < n - j - 1; ++i, ++cnt)
421 {
422 nodeList[tmp[Mode(i, j, 0)]] = tmpNodes[3][cnt];
423 nodeList[tmp[Mode(i, 0, j)]] = tmpNodes[0][cnt];
424 nodeList[tmp[Mode(n - 1 - i - j, i, j)]] = tmpNodes[1][cnt];
425 nodeList[tmp[Mode(0, i, j)]] = tmpNodes[2][cnt];
426 }
427 }
428
429 if (n == 4)
430 {
431 return nodeList;
432 }
433
434 // Finally, recurse on interior volume
435 std::vector<long long> intNodes, tmpInt;
436 for (int i = offset; i < nTet; ++i)
437 {
438 intNodes.emplace_back(nodes[i]);
439 }
440
441 tmpInt = tetTensorNodeOrdering(intNodes);
442
443 for (int k = 1, cnt = 0; k < n - 2; ++k)
444 {
445 for (int j = 1; j < n - k - 1; ++j)
446 {
447 for (int i = 1; i < n - k - j - 1; ++i)
448 {
449 nodeList[tmp[Mode(i, j, k)]] = tmpInt[cnt++];
450 }
451 }
452 }
453
454 return nodeList;
455}
456
457// This prism mapping is NOT using the standard spectral mapping! The point
458// output from ProcessEquiSpacedOutput is not the same as the rest of the code.
459std::vector<long long> prismTensorNodeOrdering(
460 const std::vector<long long> &nodes)
461{
462 int nN = static_cast<int>(nodes.size());
463 int n = static_cast<int>(std::cbrt(2 * nN));
464
465 std::vector<long long> nodeList(nN);
466 int nQuad = n * n;
467 int nPrism = n * n * (n + 1) / 2;
468 int edge = n - 2;
469
470 // Vertices
471 nodeList[0] = nodes[0];
472 if (n > 1)
473 {
474 nodeList[nPrism - n] = nodes[1];
475 nodeList[n - 1] = nodes[2];
476 nodeList[nQuad - n] = nodes[3];
477 nodeList[nPrism - 1] = nodes[4];
478 nodeList[nQuad - 1] = nodes[5];
479 }
480
481 if (n == 2)
482 {
483 return nodeList;
484 }
485
486 int nPts = 6;
487 int cnt = 0;
488 for (int i = 0; i < n - 2; ++i)
489 {
490 // Triangle 1 edges
491 cnt += n * (n - i); // This is adding the quads as they reduce a side by
492 // one each time
493 nodeList[cnt] = nodes[nPts + i]; // Edge 1
494 nodeList[cnt + edge - i] = nodes[nPts + 2 * edge - i - 1]; // Edge 2
495 nodeList[i + 1] = nodes[nPts + 3 * edge - i - 1]; // Edge 3
496 }
497
498 cnt = nQuad;
499 for (int i = 1; i < n - 1; ++i)
500 {
501 // Triangle 2 edges
502 cnt += n * (n - i); // This is adding the quads as they reduce a side by
503 // one each time
504 nodeList[cnt - n + i] = nodes[nPts + 3 * edge + i - 1]; // Edge 1 (4)
505 nodeList[cnt - 1] = nodes[nPts + 5 * edge - i]; // Edge 2 (5)
506 nodeList[nQuad - i - 1] = nodes[nPts + 5 * edge + i - 1]; // Edge 3 (6)
507 }
508
509 // Vertical edges
510 for (int i = 1; i < n - 1; ++i)
511 {
512 // Vertical prism edges
513 nodeList[n * i] = nodes[nPts + 6 * edge + i - 1]; // Edge 1 (7)
514 nodeList[nPrism - n + i] = nodes[nPts + 7 * edge + i - 1]; // Edge 2 (8)
515 nodeList[n * i + n - 1] = nodes[nPts + 8 * edge + i - 1]; // Edge 3 (9)
516 }
517
518 // Do quad faces first as when n == 3 we have no tri faces
519 int nSmallTri = (n - 3) * (n - 2) / 2;
520 int nSmallQuad = (n - 2) * (n - 2);
521 nPts = 6 + 9 * edge + 2 * nSmallTri;
522 cnt = nQuad;
523 for (int i = 1; i < n - 1; ++i)
524 {
525 for (int j = 0; j < n - 2; ++j)
526 {
527 nodeList[cnt + (n - i) + (n - i) * j] =
528 nodes[nPts + j * edge + (i - 1)];
529
530 nodeList[cnt + 2 * (n - i) - 1 + (n - i) * j] =
531 nodes[nPts + nSmallQuad + (j + 1) * edge - i];
532
533 nodeList[i * n + j + 1] =
534 nodes[nPts + 2 * nSmallQuad + i * edge - j - 1];
535 }
536
537 cnt += n * (n - i);
538 }
539
540 if (n == 3)
541 {
542 return nodeList;
543 }
544
545 // Triangular alignment
546 std::vector<long long> tmpNodes(nSmallTri);
547 std::iota(tmpNodes.begin(), tmpNodes.end(), 0);
548 if (n > 4)
549 {
550 std::vector<long long> triVertId(3), toAlign(3);
551 triVertId[0] = 0;
552 triVertId[1] = 1;
553 triVertId[2] = 2;
554
555 toAlign[0] = 0;
556 toAlign[1] = 2;
557 toAlign[2] = 1;
558 Align(triVertId, toAlign, tmpNodes);
559 }
560
561 // Triangle face 1
562 nPts = 6 + 9 * edge;
563 cnt = nQuad;
564 int cnt2 = 0;
565 for (int i = 1; i < n - 2; ++i)
566 {
567 for (int j = 1; j < n - i - 1; ++j)
568 {
569 nodeList[cnt + j] = nPts + tmpNodes[cnt2++];
570 }
571
572 cnt += n * (n - i);
573 }
574
575 // Triangle face 2
576 cnt = nQuad;
577 cnt2 = 0;
578 for (int i = 1; i < n - 2; ++i)
579 {
580 cnt += n * (n - i);
581 for (int j = 1; j < n - i - 1; ++j)
582 {
583 nodeList[cnt - (n - i) + j] = nPts + nSmallTri + tmpNodes[cnt2++];
584 }
585 }
586
587 // Triangular volume (uses same alignment mapping)
588 nPts = 6 + 9 * edge + 2 * nSmallTri + 3 * nSmallQuad;
589 for (int k = 1; k < n - 1; ++k)
590 {
591 cnt = nQuad;
592 cnt2 = 0;
593 for (int i = 1; i < n - 2; ++i)
594 {
595 for (int j = 1; j < n - i - 1; ++j)
596 {
597 nodeList[cnt + k * (n - i) + j] = nPts + tmpNodes[cnt2++];
598 }
599
600 cnt += n * (n - i);
601 }
602
603 nPts += nSmallTri;
604 }
605
606 return nodeList;
607}
608
609std::vector<long long> hexTensorNodeOrdering(
610 const std::vector<long long> &nodes)
611{
612 int nN = static_cast<int>(nodes.size());
613 int n = static_cast<int>(std::cbrt(nN));
614
615 std::vector<long long> nodeList(nN);
616
617 // Vertices
618 nodeList[0] = nodes[0];
619 if (n == 1)
620 {
621 return nodeList;
622 }
623
624 // Vertices: same order as Nektar++
625 nodeList[n - 1] = nodes[1];
626 nodeList[n * n - 1] = nodes[2];
627 nodeList[n * (n - 1)] = nodes[3];
628 nodeList[n * n * (n - 1)] = nodes[4];
629 nodeList[n - 1 + n * n * (n - 1)] = nodes[5];
630 nodeList[n * n - 1 + n * n * (n - 1)] = nodes[6];
631 nodeList[n * (n - 1) + n * n * (n - 1)] = nodes[7];
632
633 if (n == 2)
634 {
635 return nodeList;
636 }
637
638 int hexEdges[12][2] = {{0, 1},
639 {n - 1, n},
640 {n * n - 1, -1},
641 {n * (n - 1), -n},
642 {0, n * n},
643 {n - 1, n * n},
644 {n * n - 1, n * n},
645 {n * (n - 1), n * n},
646 {n * n * (n - 1), 1},
647 {n * n * (n - 1) + n - 1, n},
648 {n * n * n - 1, -1},
649 {n * n * (n - 1) + n * (n - 1), -n}};
650 int hexFaces[6][3] = {{0, 1, n}, {0, 1, n * n},
651 {n - 1, n, n * n}, {n * (n - 1), 1, n * n},
652 {0, n, n * n}, {n * n * (n - 1), 1, n}};
653
654 int gmshToNekEdge[12] = {0, 1, -2, -3, 8, 9, -10, -11, 4, 5, 6, 7};
655
656 // Edges
657 int offset = 8;
658 for (int i : gmshToNekEdge)
659 {
660 int e = abs(i);
661
662 if (i >= 0)
663 {
664 for (int j = 1; j < n - 1; ++j)
665 {
666 nodeList[hexEdges[e][0] + j * hexEdges[e][1]] = nodes[offset++];
667 }
668 }
669 else
670 {
671 for (int j = 1; j < n - 1; ++j)
672 {
673 nodeList[hexEdges[e][0] + (n - j - 1) * hexEdges[e][1]] =
674 nodes[offset++];
675 }
676 }
677 }
678
679 // Faces
680 int gmsh2NekFace[6] = {4, 2, 1, 3, 0, 5};
681
682 // Map which defines orientation between Gmsh and Nektar++ faces.
683 StdRegions::Orientation faceOrient[6] = {
690
691 for (int i = 0; i < 6; ++i)
692 {
693 int n2 = (n - 2) * (n - 2);
694 int face = gmsh2NekFace[i];
695 offset = 8 + 12 * (n - 2) + i * n2;
696
697 // Create a list of interior face nodes for this face only.
698 std::vector<long long> faceNodes(n2);
699 for (int j = 0; j < n2; ++j)
700 {
701 faceNodes[j] = nodes[offset + j];
702 }
703
704 // Now get the reordering of this face, which puts Gmsh
705 // recursive ordering into Nektar++ row-by-row order.
706 std::vector<long long> tmp(n2);
707
708 // Finally reorient the face according to the geometry
709 // differences.
710 if (faceOrient[i] == StdRegions::eDir1FwdDir1_Dir2FwdDir2)
711 {
712 // Orientation is the same, just copy.
713 tmp = faceNodes;
714 }
715 else if (faceOrient[i] == StdRegions::eDir1FwdDir2_Dir2FwdDir1)
716 {
717 // Tranposed faces
718 for (int j = 0; j < n - 2; ++j)
719 {
720 for (int k = 0; k < n - 2; ++k)
721 {
722 tmp[j * (n - 2) + k] = faceNodes[k * (n - 2) + j];
723 }
724 }
725 }
726 else if (faceOrient[i] == StdRegions::eDir1BwdDir1_Dir2FwdDir2)
727 {
728 for (int j = 0; j < n - 2; ++j)
729 {
730 for (int k = 0; k < n - 2; ++k)
731 {
732 tmp[j * (n - 2) + k] = faceNodes[j * (n - 2) + (n - k - 3)];
733 }
734 }
735 }
736
737 // Now put this into the right place in the output array
738 for (int k = 1; k < n - 1; ++k)
739 {
740 for (int j = 1; j < n - 1; ++j)
741 {
742 nodeList[hexFaces[face][0] + j * hexFaces[face][1] +
743 k * hexFaces[face][2]] =
744 faceNodes[(k - 1) * (n - 2) + j - 1];
745 }
746 }
747 }
748
749 // Finally, recurse on interior volume
750 std::vector<long long> intNodes, tmpInt;
751 for (int i = 8 + 12 * (n - 2) + 6 * (n - 2) * (n - 2); i < n * n * n; ++i)
752 {
753 intNodes.push_back(nodes[i]);
754 }
755
756 if (!intNodes.empty())
757 {
758 tmpInt = hexTensorNodeOrdering(intNodes);
759 for (int k = 1, cnt = 0; k < n - 1; ++k)
760 {
761 for (int j = 1; j < n - 1; ++j)
762 {
763 for (int i = 1; i < n - 1; ++i)
764 {
765 nodeList[i + j * n + k * n * n] = tmpInt[cnt++];
766 }
767 }
768 }
769 }
770
771 return nodeList;
772}
773
774std::vector<long long> lowOrderMapping(int nDim, Array<OneD, int> nquad)
775{
776 std::vector<long long> p;
777 switch (nDim)
778 {
779 case 3:
780 for (int j = 0; j < nquad[0] - 1; ++j)
781 {
782 for (int k = 0; k < nquad[1] - 1; ++k)
783 {
784 for (int l = 0; l < nquad[2] - 1; ++l)
785 {
786 p.insert(
787 p.end(),
788 {l * nquad[0] * nquad[1] + k * nquad[0] + j,
789 l * nquad[0] * nquad[1] + k * nquad[0] + j + 1,
790 l * nquad[0] * nquad[1] + (k + 1) * nquad[0] + j +
791 1,
792 l * nquad[0] * nquad[1] + (k + 1) * nquad[0] + j,
793 (l + 1) * nquad[0] * nquad[1] + k * nquad[0] + j,
794 (l + 1) * nquad[0] * nquad[1] + k * nquad[0] + j +
795 1,
796 (l + 1) * nquad[0] * nquad[1] +
797 (k + 1) * nquad[0] + j + 1,
798 (l + 1) * nquad[0] * nquad[1] +
799 (k + 1) * nquad[0] + j});
800 }
801 }
802 }
803 break;
804 case 2:
805 for (int j = 0; j < nquad[0] - 1; ++j)
806 {
807 for (int k = 0; k < nquad[1] - 1; ++k)
808 {
809 p.insert(p.end(), {k * nquad[0] + j, k * nquad[0] + j + 1,
810 (k + 1) * nquad[0] + j + 1,
811 (k + 1) * nquad[0] + j});
812 }
813 }
814 break;
815 case 1:
816 for (int j = 0; j < nquad[0] - 1; ++j)
817 {
818 p.insert(p.end(), {j, j + 1});
819 }
820 break;
821 default:
822 break;
823 }
824
825 return p;
826}
827
828int GetHighOrderVtkCellType([[maybe_unused]] int sType)
829{
830#if VTK_MAJOR_VERSION >= 8
831 // For deformed elements this is a map of shape type to VTK type
832 static const std::map<int, int> vtkCellType = {
833 {{LibUtilities::eSegment, VTK_LAGRANGE_CURVE},
834 {LibUtilities::eTriangle, VTK_LAGRANGE_TRIANGLE},
835 {LibUtilities::eQuadrilateral, VTK_LAGRANGE_QUADRILATERAL},
836 {LibUtilities::eTetrahedron, VTK_LAGRANGE_TETRAHEDRON},
837 {LibUtilities::ePyramid, VTK_LAGRANGE_PYRAMID},
838 {LibUtilities::ePrism, VTK_LAGRANGE_WEDGE},
839 {LibUtilities::eHexahedron, VTK_LAGRANGE_HEXAHEDRON}}};
840
841 return vtkCellType.at(sType);
842#else
843 NEKERROR(
845 "High-order VTK output requires minimum VTK library version of 8.0")
846 return 0;
847#endif
848}
849
850} // namespace
851
852vtkSmartPointer<vtkUnstructuredGrid> OutputVtk::OutputFromExpLowOrder()
853{
854 // Mesh file to output
855 vtkSmartPointer<vtkUnstructuredGrid> vtkMesh =
856 vtkSmartPointer<vtkUnstructuredGrid>::New();
857
858 // Cache ordering so we aren't recreating mappings
859 std::unordered_map<size_t, std::vector<long long>> mappingCache;
860
861 // Choose cell types and num quad points based on dimension
862 int nDim = m_f->m_exp[0]->GetShapeDimension();
863 int nHomoDir = m_f->m_numHomogeneousDir;
864
865 auto type = (nDim + nHomoDir == 3) ? VTK_HEXAHEDRON
866 : (nDim + nHomoDir == 2) ? VTK_QUAD
867 : VTK_LINE;
868 int nQpts = (nDim + nHomoDir == 3) ? 8 : (nDim + nHomoDir == 2) ? 4 : 2;
869
870 // Fill homogeneous info
871 if (nHomoDir != 0)
872 {
873 ASSERTL0(nHomoDir == 1 &&
874 m_f->m_exp[0]->GetExpType() == MultiRegions::e3DH1D,
875 "Only regular expansions and the 3DH1D homogeneous expansion "
876 "are supported in the new VTK writer. Please use the 'legacy' "
877 "option for all other expansion types.")
878
879 m_numPlanes = m_f->m_exp[0]->GetHomogeneousBasis()->GetNumModes();
880 // Extra plane if fourier
881 m_extraPlane = (m_f->m_exp[0]->GetHomogeneousBasis()->GetBasisType() ==
883 m_f->m_exp[0]->GetHomogeneousBasis()->GetPointsType() ==
885 }
886
887 // Insert points
888 vtkNew<vtkPoints> vtkPoints;
889 vtkPoints->SetDataType(VTK_DOUBLE);
890
891 int offset = 0;
892 int nElmts = static_cast<int>(m_f->m_exp[0]->GetNumElmts());
893 for (int i = 0; i < nElmts; ++i)
894 {
895 Array<OneD, int> nquad(3, 0);
896 for (int j = 0; j < nDim; ++j)
897 {
898 nquad[j] = m_f->m_exp[0]->GetExp(i)->GetNumPoints(j);
899 }
900
901 // Add homogeneous direction num points
902 if (nHomoDir != 0)
903 {
904 nquad[nDim] = m_numPlanes;
905
906 // Add extra plane if fourier
907 if (m_extraPlane)
908 {
909 nquad[nDim]++;
910 }
911 }
912
913 // Get total points by multiplying all nquads together (accounts for
914 // homo)
915 int nPts = nquad[0];
916 for (int j = 1; j < nDim + nHomoDir; ++j)
917 {
918 nPts *= nquad[j];
919 }
920
921 // Get all coordinates at quadrature points in this element
922 Array<OneD, NekDouble> x(nPts, 0.0), y(nPts, 0.0), z(nPts, 0.0);
923 m_f->m_exp[0]->GetCoords(i, x, y, z);
924
925 // If add extra plane for fourier need to fill last plane coordinates
926 if (m_extraPlane)
927 {
928 // Copy x & y to extra plane
930 Vmath::Vcopy(nquad[0] * nquad[1], x, 1,
931 tmp = x + (nquad[2] - 1) * nquad[0] * nquad[1], 1);
932 Vmath::Vcopy(nquad[0] * nquad[1], y, 1,
933 tmp = y + (nquad[2] - 1) * nquad[0] * nquad[1], 1);
934
935 // Fill z on extra plane
936 NekDouble zDouble = z[nquad[0] * nquad[1] * (nquad[2] - 1) - 1] +
937 (z[nquad[0] * nquad[1]] - z[0]);
938 Vmath::Fill(nquad[0] * nquad[1], zDouble,
939 tmp = z + (nquad[2] - 1) * nquad[0] * nquad[1], 1);
940 }
941
942 // Insert points in to vtk mesh
943 for (int j = 0; j < nPts; ++j)
944 {
945 vtkPoints->InsertNextPoint(x[j], y[j], z[j]);
946 }
947
948 // Insert elements
949 auto oIt = mappingCache.find(key3(nquad[0], nquad[1], nquad[2]));
950 if (oIt == mappingCache.end())
951 {
952 auto p = lowOrderMapping(nDim + nHomoDir, nquad);
953 oIt = mappingCache
954 .insert(
955 std::make_pair(key3(nquad[0], nquad[1], nquad[2]), p))
956 .first;
957 }
958
959 auto p = oIt->second;
960 std::for_each(p.begin(), p.end(),
961 [offset](long long &d) { d += offset; });
962 for (int j = 0; j < p.size(); j += nQpts)
963 {
964 vtkMesh->InsertNextCell(type, nQpts, &p[j]);
965 }
966
967 offset += nPts;
968 }
969
970 vtkMesh->SetPoints(vtkPoints.GetPointer());
971
972 return vtkMesh;
973}
974
976 po::variables_map &vm, std::string &filename,
977 vtkSmartPointer<vtkUnstructuredGrid> &vtkMesh)
978{
979 // Insert field information we iterate by element from first plane to last,
980 // while the fld file iterates by plane from first element to last
981 int nElmts = static_cast<int>(m_f->m_exp[0]->GetNumElmts());
982 int numPtsPerPlane = m_f->m_exp[0]->GetTotPoints() / m_numPlanes;
983 for (int v = 0; v < m_f->m_variables.size(); ++v)
984 {
985 vtkNew<vtkDoubleArray> fieldData;
986 for (int i = 0; i < nElmts; ++i)
987 {
988 int elmtOffset = m_f->m_exp[v]->GetPhys_Offset(i);
989 int nPtsPerElmtPerPlane = m_f->m_exp[v]->GetExp(i)->GetTotPoints();
990
991 for (int j = 0; j < m_numPlanes; ++j)
992 {
993 int planeOffset = j * numPtsPerPlane;
994 for (int k = 0; k < nPtsPerElmtPerPlane; ++k)
995 {
996 fieldData->InsertNextValue(
997 m_f->m_exp[v]->GetPhys()[elmtOffset + planeOffset + k]);
998 }
999 }
1000
1001 // if extra plane we copy the first plane values in to last plane
1002 if (m_extraPlane)
1003 {
1004 for (int k = 0; k < nPtsPerElmtPerPlane; ++k)
1005 {
1006 fieldData->InsertNextValue(
1007 m_f->m_exp[v]->GetPhys()[elmtOffset + k]);
1008 }
1009 }
1010 }
1011
1012 fieldData->SetName(&m_f->m_variables[v][0]);
1013 vtkMesh->GetPointData()->AddArray(fieldData.GetPointer());
1014 }
1015
1016 WriteVTK(vtkMesh, filename, vm);
1017}
1018
1020 std::string &filename)
1021{
1022 ASSERTL0(
1023 m_f->m_numHomogeneousDir == 0,
1024 "Multi block VTK is not implemented for homogeneous expansion types.")
1025
1026 ASSERTL0(m_f->m_comm->GetSpaceComm()->IsSerial(),
1027 "Multi block VTK is not implemented in parallel.")
1028
1029 int dim = m_f->m_graph->GetMeshDimension();
1030
1031 // Create mappings from geometry id to expansion ids
1032 std::array<std::map<int, std::pair<int, int>>, 4> geomIdToExpId;
1033 int nElmts = static_cast<int>(m_f->m_exp[0]->GetNumElmts());
1034 for (int i = 0; i < nElmts; ++i)
1035 {
1036 auto geom = m_f->m_exp[0]->GetExp(i)->GetGeom();
1037 geomIdToExpId[geom->GetShapeDim()][geom->GetGlobalID()] =
1038 std::make_pair(i, -1);
1039
1040 for (int j = 0; j < geom->GetNumFaces(); ++j)
1041 {
1042 geomIdToExpId[2][geom->GetFid(j)] = std::make_pair(i, j);
1043 }
1044
1045 for (int j = 0; j < geom->GetNumEdges(); ++j)
1046 {
1047 geomIdToExpId[1][geom->GetEid(j)] = std::make_pair(i, j);
1048 }
1049 }
1050
1051 // Cache ordering so we aren't recreating mappings
1052 std::unordered_map<size_t, std::vector<long long>> mappingCache;
1053
1054 std::map<int, vtkNew<vtkUnstructuredGrid>> vtkMesh;
1055 std::map<int, SpatialDomains::CompositeSharedPtr> composites =
1056 m_f->m_graph->GetComposites();
1057 std::map<int, std::string> compositeNames;
1058 for (auto &comp : composites)
1059 {
1060 // Vector of field data
1061 std::vector<vtkNew<vtkDoubleArray>> fieldData(m_f->m_variables.size());
1062
1063 // Insert points
1064 vtkNew<vtkPoints> vtkPoints;
1065 vtkPoints->SetDataType(VTK_DOUBLE);
1066
1067 int compId = comp.first;
1068 std::vector<std::shared_ptr<SpatialDomains::Geometry>> geomVec =
1069 comp.second->m_geomVec;
1070
1071 unsigned int offset = 0;
1072 for (auto &geom : geomVec)
1073 {
1074 int geomId = geom->GetGlobalID();
1075 int sDim = geom->GetShapeDim();
1076 auto type = (sDim == 3) ? VTK_HEXAHEDRON
1077 : (sDim == 2) ? VTK_QUAD
1078 : VTK_LINE;
1079 int nQpts = (sDim == 3) ? 8 : (sDim == 2) ? 4 : 2;
1080
1082 m_f->m_exp[0]->GetExp(geomIdToExpId[sDim][geomId].first);
1083
1084 unsigned int nPts = exp->GetTotPoints();
1085 Array<OneD, NekDouble> x(nPts, 0.0), y(nPts, 0.0), z(nPts, 0.0);
1086 exp->GetCoords(x, y, z);
1087
1088 int offsetPhys = m_f->m_exp[0]->GetPhys_Offset(
1089 geomIdToExpId[sDim][geomId].first);
1090 if (sDim == dim)
1091 {
1092 for (int j = 0; j < nPts; ++j)
1093 {
1094 vtkPoints->InsertNextPoint(x[j], y[j], z[j]);
1095
1096 // Add field data
1097 for (int k = 0; k < m_f->m_variables.size(); ++k)
1098 {
1099 fieldData[k]->InsertNextValue(
1100 m_f->m_exp[k]->GetPhys()[j + offsetPhys]);
1101 }
1102 }
1103 }
1104 else
1105 {
1106 Array<OneD, int> pointsMap;
1107 exp->GetTracePhysMap(geomIdToExpId[sDim][geomId].second,
1108 pointsMap);
1109 for (int j : pointsMap)
1110 {
1111 vtkPoints->InsertNextPoint(x[j], y[j], z[j]);
1112
1113 // Add field data
1114 for (int k = 0; k < m_f->m_variables.size(); ++k)
1115 {
1116 fieldData[k]->InsertNextValue(
1117 m_f->m_exp[k]->GetPhys()[offsetPhys + j]);
1118 }
1119 }
1120
1121 exp = exp->GetTraceExp(geomIdToExpId[sDim][geomId].second);
1122 nPts = pointsMap.size();
1123 }
1124
1125 Array<OneD, int> nquad(3, 0);
1126 for (int j = 0; j < sDim; ++j)
1127 {
1128 nquad[j] = exp->GetNumPoints(j);
1129 }
1130
1131 auto oIt = mappingCache.find(key3(nquad[0], nquad[1], nquad[2]));
1132 if (oIt == mappingCache.end())
1133 {
1134 auto p = lowOrderMapping(sDim, nquad);
1135 oIt = mappingCache
1136 .insert(std::make_pair(
1137 key3(nquad[0], nquad[1], nquad[2]), p))
1138 .first;
1139 }
1140
1141 auto p = oIt->second;
1142 std::for_each(p.begin(), p.end(),
1143 [offset](long long &d) { d += offset; });
1144
1145 for (int j = 0; j < p.size(); j += nQpts)
1146 {
1147 vtkMesh[compId]->InsertNextCell(type, nQpts, &p[j]);
1148 }
1149
1150 offset += nPts;
1151 }
1152
1153 vtkMesh[compId]->SetPoints(vtkPoints.GetPointer());
1154
1155 // Add all fields to vtkMesh
1156 for (int k = 0; k < m_f->m_variables.size(); ++k)
1157 {
1158 fieldData[k]->SetName(&m_f->m_variables[k][0]);
1159 vtkMesh[compId]->GetPointData()->AddArray(
1160 fieldData[k].GetPointer());
1161 }
1162
1163 // Find composite names if exists and store in map
1164 // Set name as boundary label if it exists otherwise index ID
1165 compositeNames[compId] = "Composite ID " + std::to_string(compId);
1166 auto clabels = m_f->m_graph->GetCompositesLabels();
1167 auto oIt = clabels.find(compId);
1168 if (oIt != clabels.end())
1169 {
1170 compositeNames[compId] = oIt->second;
1171 }
1172 }
1173
1174 // Main multi-block set
1175 vtkNew<vtkMultiBlockDataSet> mainBlock;
1176
1177 std::set<int> compSet;
1178
1179 // Fill domain multi blocks from composites
1180 vtkNew<vtkMultiBlockDataSet> mainDomainBlock;
1181 auto domains = m_f->m_graph->GetDomain();
1182 std::vector<vtkNew<vtkMultiBlockDataSet>> domainMultiBlocks(domains.size());
1183 for (int i = 0; i < domains.size(); ++i)
1184 {
1185 auto dom = domains[i];
1186
1187 // Loop over composites and see if in domain
1188 for (auto &comp : composites)
1189 {
1190 int compId = comp.first;
1191 if (dom.find(compId) != dom.end())
1192 {
1193 unsigned int nBlock = domainMultiBlocks[i]->GetNumberOfBlocks();
1194 domainMultiBlocks[i]->SetBlock(nBlock,
1195 vtkMesh[compId].GetPointer());
1196 domainMultiBlocks[i]->GetMetaData(nBlock)->Set(
1197 vtkCompositeDataSet::NAME(),
1198 compositeNames[compId].c_str());
1199 compSet.insert(compId);
1200 }
1201 }
1202
1203 unsigned int nBlock = mainDomainBlock->GetNumberOfBlocks();
1204 mainDomainBlock->SetBlock(nBlock, domainMultiBlocks[i].GetPointer());
1205 mainDomainBlock->GetMetaData(nBlock)->Set(
1206 vtkCompositeDataSet::NAME(),
1207 ("Domain ID " + std::to_string(i)).c_str());
1208 }
1209
1210 if (mainDomainBlock->GetNumberOfBlocks() != 0)
1211 {
1212 auto nBlock = static_cast<unsigned int>(mainBlock->GetNumberOfBlocks());
1213 mainBlock->SetBlock(nBlock, mainDomainBlock.GetPointer());
1214 mainBlock->GetMetaData(nBlock)->Set(vtkCompositeDataSet::NAME(),
1215 "Domains");
1216 }
1217
1218 // Fill boundary multi blocks from composites
1220 m_f->m_exp[0]->GetGraph());
1222 bcs.GetBoundaryRegions();
1223
1224 vtkNew<vtkMultiBlockDataSet> mainBoundaryBlock;
1225 std::vector<vtkNew<vtkMultiBlockDataSet>> boundaryMultiBlocks(
1226 bregions.size());
1227 int cnt = 0;
1228 for (auto &boundary : bregions)
1229 {
1230 // Loop over composites and see if in boundary
1231 for (auto &comp : composites)
1232 {
1233 int compId = comp.first;
1234 if (boundary.second->find(compId) != boundary.second->end())
1235 {
1236 unsigned int nBlock =
1237 boundaryMultiBlocks[cnt]->GetNumberOfBlocks();
1238 boundaryMultiBlocks[cnt]->SetBlock(
1239 nBlock, vtkMesh[compId].GetPointer());
1240 boundaryMultiBlocks[cnt]->GetMetaData(nBlock)->Set(
1241 vtkCompositeDataSet::NAME(),
1242 compositeNames[compId].c_str());
1243 compSet.insert(compId);
1244 }
1245 }
1246
1247 unsigned int nBlock = mainBoundaryBlock->GetNumberOfBlocks();
1248 mainBoundaryBlock->SetBlock(nBlock,
1249 boundaryMultiBlocks[cnt++].GetPointer());
1250
1251 // Set name as boundary label if it exists otherwise index ID
1252 std::string name = "Boundary ID " + std::to_string(boundary.first);
1253 auto blabels = bcs.GetBoundaryLabels();
1254 auto oIt = blabels.find(boundary.first);
1255 if (oIt != blabels.end())
1256 {
1257 name = oIt->second;
1258 }
1259
1260 mainBoundaryBlock->GetMetaData(nBlock)->Set(vtkCompositeDataSet::NAME(),
1261 name.c_str());
1262 }
1263
1264 if (mainBoundaryBlock->GetNumberOfBlocks() != 0)
1265 {
1266 auto nBlock = static_cast<unsigned int>(mainBlock->GetNumberOfBlocks());
1267 mainBlock->SetBlock(nBlock, mainBoundaryBlock.GetPointer());
1268 mainBlock->GetMetaData(nBlock)->Set(vtkCompositeDataSet::NAME(),
1269 "Boundaries");
1270 }
1271
1272 // Include all other composites (not domains or boundaries)
1273 vtkNew<vtkMultiBlockDataSet> compsMultiBlocks;
1274 cnt = 0;
1275 for (auto &comp : composites)
1276 {
1277 int compId = comp.first;
1278 if (compSet.find(compId) == compSet.end())
1279 {
1280 unsigned int nBlock = compsMultiBlocks->GetNumberOfBlocks();
1281 compsMultiBlocks->SetBlock(nBlock, vtkMesh[compId].GetPointer());
1282 compsMultiBlocks->GetMetaData(nBlock)->Set(
1283 vtkCompositeDataSet::NAME(), compositeNames[compId].c_str());
1284 }
1285 }
1286
1287 if (compsMultiBlocks->GetNumberOfBlocks() != 0)
1288 {
1289 auto nBlock = static_cast<unsigned int>(mainBlock->GetNumberOfBlocks());
1290 mainBlock->SetBlock(nBlock, compsMultiBlocks.GetPointer());
1291 mainBlock->GetMetaData(nBlock)->Set(vtkCompositeDataSet::NAME(),
1292 "Other composites");
1293 }
1294
1295 WriteVTK(mainBlock.GetPointer(), filename, vm);
1296}
1297
1298vtkSmartPointer<vtkUnstructuredGrid> OutputVtk::OutputFromExpHighOrder(
1299 po::variables_map &vm)
1300{
1301 ASSERTL0(
1302 m_f->m_numHomogeneousDir == 0,
1303 "High order VTK is not implemented for homogeneous expansion types.")
1304
1305 // Save shapetypes before convert to equispaced because that nukes the
1306 // explist
1307 std::vector<LibUtilities::ShapeType> sType;
1308 for (auto &i : *m_f->m_exp[0]->GetExp())
1309 {
1310 sType.emplace_back(i->GetGeom()->GetShapeType());
1311 }
1312
1313 // Convert expansion to an equispaced points field
1314 auto equispaced =
1316 equispaced->Process(vm);
1317
1318 // Insert points
1319 vtkNew<vtkPoints> vtkPoints;
1320 LibUtilities::PtsFieldSharedPtr fPts = m_f->m_fieldPts;
1321 vtkPoints->SetDataType(VTK_DOUBLE);
1322
1324 fPts->GetPts(pts);
1325
1326 int dim = static_cast<int>(fPts->GetDim());
1327 int nPts = static_cast<int>(fPts->GetNpoints());
1328 switch (dim)
1329 {
1330 case 3:
1331 for (int i = 0; i < nPts; ++i)
1332 {
1333 vtkPoints->InsertNextPoint(pts[0][i], pts[1][i], pts[2][i]);
1334 }
1335 break;
1336 case 2:
1337 for (int i = 0; i < nPts; ++i)
1338 {
1339 vtkPoints->InsertNextPoint(pts[0][i], pts[1][i], 0.0);
1340 }
1341 break;
1342 case 1:
1343 for (int i = 0; i < nPts; ++i)
1344 {
1345 vtkPoints->InsertNextPoint(pts[0][i], 0.0, 0.0);
1346 }
1347 break;
1348 default:
1349 break;
1350 }
1351
1352 // Mesh file to output
1353 vtkSmartPointer<vtkUnstructuredGrid> vtkMesh =
1354 vtkSmartPointer<vtkUnstructuredGrid>::New();
1355 vtkMesh->SetPoints(vtkPoints.GetPointer());
1356
1357 // Cache ordering for shape type & npts so we aren't recreating mappings
1358 std::unordered_map<size_t, std::vector<long long>> mappingCache;
1359
1360 // Get offset per element in a vector from ppe
1361 std::vector<int> ppe = m_f->m_fieldPts->GetPointsPerElement();
1362 Array<OneD, int> ppeOffset(ppe.size() + 1, 0.0);
1363 std::partial_sum(ppe.begin(), ppe.end(), &ppeOffset[1]);
1364
1365 // Insert elements
1366 for (int i = 0; i < ppe.size(); ++i)
1367 {
1368 // Construct inverse of input reordering.
1369 // First try to find it in our mapping cache.
1370 auto oIt = mappingCache.find(key2(sType[i], ppe[i]));
1371 if (oIt == mappingCache.end())
1372 {
1373 std::vector<long long> p(ppe[i]);
1374 std::iota(p.begin(), p.end(), 0);
1375 switch (sType[i])
1376 {
1378 p = quadTensorNodeOrdering(p);
1379 break;
1381 p = triTensorNodeOrdering(p);
1382 break;
1384 p = tetTensorNodeOrdering(p);
1385 break;
1387 p = hexTensorNodeOrdering(p);
1388 break;
1390 p = prismTensorNodeOrdering(p);
1391 break;
1392 default:
1394 "High-order VTU output not set up for the " +
1395 static_cast<std::string>(
1396 LibUtilities::ShapeTypeMap[sType[i]]) +
1397 " shape type.");
1398 break;
1399 }
1400
1401 // Invert the ordering as this is for spectral -> recursive (VTU)
1402 std::vector<long long> inv(ppe[i]);
1403 for (int j = 0; j < ppe[i]; ++j)
1404 {
1405 inv[p[j]] = j;
1406 }
1407
1408 oIt =
1409 mappingCache.insert(std::make_pair(key2(sType[i], ppe[i]), inv))
1410 .first;
1411 }
1412
1413 // Add offset to reordering
1414 std::vector<long long> p = oIt->second;
1415 for (long long &j : p)
1416 {
1417 j += ppeOffset[i];
1418 }
1419
1420 vtkMesh->InsertNextCell(GetHighOrderVtkCellType(sType[i]), ppe[i],
1421 &p[0]);
1422 }
1423
1424 return vtkMesh;
1425}
1426
1428 po::variables_map &vm, std::string &filename,
1429 vtkSmartPointer<vtkUnstructuredGrid> &vtkMesh)
1430{
1431 // Convert expansion to an equispaced points field
1432 if (m_cachedMesh)
1433 {
1434 auto equispaced =
1436 equispaced->Process(vm);
1437 }
1438 LibUtilities::PtsFieldSharedPtr fPts = m_f->m_fieldPts;
1439
1440 int nPts = static_cast<int>(fPts->GetNpoints());
1441 int dim = static_cast<int>(fPts->GetDim());
1442
1444 fPts->GetPts(pts);
1445
1446 // Insert field information
1447 for (int i = 0; i < fPts->GetNFields(); ++i)
1448 {
1449 vtkNew<vtkDoubleArray> fieldData;
1450 fieldData->SetArray(&pts[dim + i][0], nPts, 1);
1451 fieldData->SetName(&fPts->GetFieldName(i)[0]);
1452 vtkMesh->GetPointData()->AddArray(fieldData.GetPointer());
1453 }
1454
1455 WriteVTK(vtkMesh, filename, vm);
1456}
1457
1458void OutputVtk::WriteVTK(vtkDataObject *vtkMesh, std::string &filename,
1459 po::variables_map &vm)
1460{
1461 // Initialise base writer class as default structured grid writer
1462 vtkSmartPointer<vtkXMLWriter> vtkMeshWriter =
1463 vtkNew<vtkXMLUnstructuredGridWriter>().GetPointer();
1464
1465 // If multi block we need to switch to the MultiBlock writer
1466 if (m_config["multiblock"].m_beenSet)
1467 {
1468 vtkMeshWriter = vtkNew<vtkXMLMultiBlockDataWriter>().GetPointer();
1469 }
1470
1471 // Choose the correct extension based on which mesh writer is used
1472 int dot = static_cast<int>(filename.find_last_of('.'));
1473 std::string body = filename.substr(0, dot + 1);
1474 filename = body + vtkMeshWriter->GetDefaultFileExtension();
1475
1476 // Set time information if exists
1477 if (!m_f->m_fieldMetaDataMap["Time"].empty())
1478 {
1479 vtkMesh->GetInformation()->Set(
1480 vtkDataObject::DATA_TIME_STEP(),
1481 std::stod(m_f->m_fieldMetaDataMap["Time"]));
1482 }
1483 else if (!m_f->m_fieldMetaDataMap["FinalTime"].empty())
1484 {
1485 vtkMesh->GetInformation()->Set(
1486 vtkDataObject::DATA_TIME_STEP(),
1487 std::stod(m_f->m_fieldMetaDataMap["FinalTime"]));
1488 }
1489
1490 // Write out the new mesh in XML format (don't support legacy
1491 // format here as we still have standard OutputVtk.cpp)
1492 vtkMeshWriter->SetFileName(filename.c_str());
1493
1494#if VTK_MAJOR_VERSION <= 5
1495 vtkMeshWriter->SetInput(vtkMesh);
1496#else
1497 vtkMeshWriter->SetInputData(vtkMesh);
1498#endif
1499
1500 if (m_config["uncompress"].m_beenSet)
1501 {
1502 vtkMeshWriter->SetDataModeToAscii();
1503 }
1504
1505 vtkMeshWriter->Update();
1506
1507 std::cout << "Written file: " << filename << std::endl;
1508
1509 // Output parallel file information if needed - using a manual write.
1510 // We could use the VTK lib to do this, but that requires VTK with MPI
1511 // enabled & messing about with the parallel controller & changing
1512 // our file naming scheme as VTK forces _${proc-number} as a suffix...
1513 if (m_f->m_comm->GetSpaceComm()->TreatAsRankZero() &&
1514 !m_f->m_comm->GetSpaceComm()->IsSerial())
1515 {
1516 WritePVtu(vm);
1517 }
1518}
1519
1520void OutputVtk::WritePVtu(po::variables_map &vm)
1521{
1522 std::string filename = m_config["outfile"].as<std::string>();
1523 int dot = static_cast<int>(filename.find_last_of('.'));
1524 std::string body = filename.substr(0, dot);
1525 filename = body + ".pvtu";
1526
1527 std::ofstream outfile(filename.c_str());
1528
1529 int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
1530 std::string path =
1532
1533 outfile << "<?xml version=\"1.0\"?>" << endl;
1534 outfile << R"(<VTKFile type="PUnstructuredGrid" version="0.1" )"
1535 << "byte_order=\"LittleEndian\">" << endl;
1536 outfile << " <PUnstructuredGrid GhostLevel=\"0\">" << endl;
1537
1538 // Add time if exists
1539 if (!m_f->m_fieldMetaDataMap["Time"].empty())
1540 {
1541 outfile << " <FieldData> " << endl;
1542 outfile << " <DataArray type=\"Float64\" Name=\"TimeValue\" "
1543 "NumberOfTuples=\"1\" format=\"ascii\">"
1544 << endl;
1545 outfile << m_f->m_fieldMetaDataMap["Time"] << std::endl;
1546 outfile << " </DataArray>" << std::endl;
1547 outfile << " </FieldData> " << endl;
1548 }
1549
1550 // Add point coordinates
1551 outfile << " <PPoints> " << endl;
1552 outfile << R"( <PDataArray type="Float64" NumberOfComponents=")" << 3
1553 << "\"/> " << endl;
1554 outfile << " </PPoints>" << endl;
1555
1556 // Add cell information
1557 outfile << " <PCells>" << endl;
1558 outfile << " <PDataArray type=\"Int64\" Name=\"connectivity\" "
1559 "NumberOfComponents=\"1\"/>"
1560 << endl;
1561 outfile << " <PDataArray type=\"Int64\" Name=\"offsets\" "
1562 "NumberOfComponents=\"1\"/>"
1563 << endl;
1564 outfile << " <PDataArray type=\"UInt8\" Name=\"types\" "
1565 "NumberOfComponents=\"1\"/>"
1566 << endl;
1567 outfile << " </PCells>" << endl;
1568
1569 // Add fields (point data)
1570 outfile << " <PPointData>" << endl;
1571 for (auto &var : m_f->m_variables)
1572 {
1573 outfile << R"( <PDataArray type="Float64" Name=")" << var << "\"/>"
1574 << endl;
1575 }
1576 outfile << " </PPointData>" << endl;
1577
1578 // Add parallel files
1579 for (int i = 0; i < nprocs; ++i)
1580 {
1581 boost::format pad("P%1$07d.vtu");
1582 pad % i;
1583 outfile << " <Piece Source=\"" << path << "/" << pad.str() << "\"/>"
1584 << endl;
1585 }
1586 outfile << " </PUnstructuredGrid>" << endl;
1587 outfile << "</VTKFile>" << endl;
1588
1589 cout << "Written file: " << filename << endl;
1590}
1591
1592void OutputVtk::v_OutputFromData([[maybe_unused]] po::variables_map &vm)
1593{
1595 "OutputVtk can't write using only FieldData. You may need "
1596 "to add a mesh XML file to your input files.");
1597}
1598
1599void OutputVtk::v_OutputFromPts(po::variables_map &vm)
1600{
1602}
1603
1604void OutputVtk::v_OutputFromExp(po::variables_map &vm)
1605{
1606 if (m_config["legacy"].m_beenSet)
1607 {
1608 ASSERTL0(!m_config["multiblock"].m_beenSet,
1609 "Multi block VTK is not implemented for legacy output.")
1610
1611 ASSERTL0(!m_config["highorder"].m_beenSet,
1612 "High order VTK is not implemented for legacy output.")
1613
1614 // No caching of mesh data in legacy output
1616 return;
1617 }
1618
1619 // Extract the output filename and extension
1620 std::string filename = OutputVtkBase::PrepareOutput(vm);
1621
1622 // Save mesh state (if using filter this allows us to only ProcessEquispaced
1623 // if needed)
1625
1626 if (m_config["highorder"].m_beenSet)
1627 {
1628 ASSERTL0(!m_config["multiblock"].m_beenSet,
1629 "Multi block VTK is not implemented for high-order output.")
1630
1631 if (!m_cachedMesh)
1632 {
1634 }
1635
1637 }
1638 else if (m_config["multiblock"].m_beenSet)
1639 {
1640 // No mesh caching for multiblock due to complexity of field data
1641 OutputFromExpLowOrderMultiBlock(vm, filename);
1642 }
1643 else
1644 {
1645 if (!m_cachedMesh)
1646 {
1648 }
1649
1650 AddFieldDataToVTKLowOrder(vm, filename, m_vtkMesh);
1651 }
1652}
1653
1654} // namespace Nektar::FieldUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272
Converter from fld to vtk.
Definition: OutputVtkBase.h:46
std::string PrepareOutput(po::variables_map &vm)
fs::path v_GetPath(std::string &filename, po::variables_map &vm) override
void v_OutputFromExp(po::variables_map &vm) override
Write from m_exp to output file.
void v_OutputFromPts(po::variables_map &vm) override
Write from pts to output file.
vtkSmartPointer< vtkUnstructuredGrid > OutputFromExpLowOrder()
Prepare low order VTK output.
Definition: OutputVtk.cpp:852
void AddFieldDataToVTKHighOrder(po::variables_map &vm, std::string &filename, vtkSmartPointer< vtkUnstructuredGrid > &vtkMesh)
Add field data to high order Lagrange VTK output.
Definition: OutputVtk.cpp:1427
void OutputFromExpLowOrderMultiBlock(po::variables_map &vm, std::string &filename)
Prepare low order multi-block VTK output & add field data.
Definition: OutputVtk.cpp:1019
OutputVtk(FieldSharedPtr f)
Definition: OutputVtk.cpp:56
vtkSmartPointer< vtkUnstructuredGrid > m_vtkMesh
Cache file for unstructured grid VTK mesh data.
Definition: OutputVtk.h:80
void v_OutputFromExp(po::variables_map &vm) final
Write from m_exp to output file.
Definition: OutputVtk.cpp:1604
void v_OutputFromPts(po::variables_map &vm) final
Write from pts to output file.
Definition: OutputVtk.cpp:1599
bool m_extraPlane
Flag if extra plane in case of fourier expansion in homogeneous dir.
Definition: OutputVtk.h:86
int m_numPlanes
Number of planes if homogeneous.
Definition: OutputVtk.h:83
void v_OutputFromData(po::variables_map &vm) final
Write from data to output file.
Definition: OutputVtk.cpp:1592
void AddFieldDataToVTKLowOrder(po::variables_map &vm, std::string &filename, vtkSmartPointer< vtkUnstructuredGrid > &vtkMesh)
Add field data to low order VTK output.
Definition: OutputVtk.cpp:975
void WriteVTK(vtkDataObject *vtkMesh, std::string &filename, po::variables_map &vm)
Write VTK file using.
Definition: OutputVtk.cpp:1458
vtkSmartPointer< vtkUnstructuredGrid > OutputFromExpHighOrder(po::variables_map &vm)
Prepare high order Lagrange VTK output.
Definition: OutputVtk.cpp:1298
static ModuleKey m_className
Definition: OutputVtk.h:58
bool m_cachedMesh
Flag if mesh has been cached.
Definition: OutputVtk.h:89
static std::shared_ptr< Module > create(const FieldSharedPtr &f)
Creates an instance of this class.
Definition: OutputVtk.h:53
void WritePVtu(po::variables_map &vm)
Write the parallel .pvtu file.
Definition: OutputVtk.cpp:1520
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
const std::map< int, std::string > & GetBoundaryLabels(void) const
Definition: Conditions.h:265
const BoundaryRegionCollection & GetBoundaryRegions(void) const
Definition: Conditions.h:234
def copy(self)
Definition: pycml.py:2663
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:990
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:75
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
Definition: Filesystem.hpp:56
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:184
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:210
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
double NekDouble
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294
Represents a command-line configuration option.
Definition: Module.h:129