Nektar++
OutputVtk.h
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: OutputVtk.h
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 output module using the VTK library
32//
33////////////////////////////////////////////////////////////////////////////////
34
35#ifndef FIELDUTILS_OUTPUTVTK
36#define FIELDUTILS_OUTPUTVTK
37
38#include "OutputVtkBase.h"
39#include <tinyxml.h>
40
41#include <vtkNew.h>
42#include <vtkSmartPointer.h>
43#include <vtkUnstructuredGrid.h>
44
45namespace Nektar::FieldUtils
46{
47
48/// Converter from fld to vtk.
49class OutputVtk final : public OutputVtkBase
50{
51public:
52 /// Creates an instance of this class
53 static std::shared_ptr<Module> create(const FieldSharedPtr &f)
54 {
56 }
57
59
60 explicit OutputVtk(FieldSharedPtr f);
61
62 ~OutputVtk() final = default;
63
64 vtkUnstructuredGrid *GetVtkGrid()
65 {
66 return m_vtkMesh.GetPointer();
67 }
68
69protected:
70 std::string v_GetModuleName() final
71 {
72 return "OutputVtk";
73 }
74
75 /// Write from pts to output file.
76 void v_OutputFromPts(po::variables_map &vm) final;
77
78 /// Write from m_exp to output file.
79 void v_OutputFromExp(po::variables_map &vm) final;
80
81 /// Write from data to output file.
82 void v_OutputFromData(po::variables_map &vm) final;
83
84 /// Cache file for unstructured grid VTK mesh data
85 vtkSmartPointer<vtkUnstructuredGrid> m_vtkMesh;
86
87 /// Number of planes if homogeneous
88 int m_numPlanes = 1;
89
90 /// Flag if extra plane in case of fourier expansion in homogeneous dir
91 bool m_extraPlane = false;
92
93 /// Flag if mesh has been cached
94 bool m_cachedMesh = false;
95
96private:
97 /// Prepare high order Lagrange VTK output
98 vtkSmartPointer<vtkUnstructuredGrid> OutputFromExpHighOrder(
99 po::variables_map &vm);
100
101 /// Add field data to high order Lagrange VTK output
103 po::variables_map &vm, std::string &filename,
104 vtkSmartPointer<vtkUnstructuredGrid> &vtkMesh);
105
106 /// Prepare low order VTK output
107 vtkSmartPointer<vtkUnstructuredGrid> OutputFromExpLowOrder();
108
109 /// Add field data to low order VTK output
111 po::variables_map &vm, std::string &filename,
112 vtkSmartPointer<vtkUnstructuredGrid> &vtkMesh);
113
114 /// Prepare low order multi-block VTK output & add field data
115 void OutputFromExpLowOrderMultiBlock(po::variables_map &vm,
116 std::string &filename);
117
118 /// Write VTK file using @param vtkMesh
119 void WriteVTK(vtkDataObject *vtkMesh, std::string &filename,
120 po::variables_map &vm);
121
122 /// Write the parallel .pvtu file
123 void WritePVtu(po::variables_map &vm);
124};
125} // namespace Nektar::FieldUtils
126
127#endif
Converter from fld to vtk.
Definition: OutputVtkBase.h:46
Converter from fld to vtk.
Definition: OutputVtk.h:50
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:1433
void OutputFromExpLowOrderMultiBlock(po::variables_map &vm, std::string &filename)
Prepare low order multi-block VTK output & add field data.
Definition: OutputVtk.cpp:1022
OutputVtk(FieldSharedPtr f)
Definition: OutputVtk.cpp:56
vtkSmartPointer< vtkUnstructuredGrid > m_vtkMesh
Cache file for unstructured grid VTK mesh data.
Definition: OutputVtk.h:85
void v_OutputFromExp(po::variables_map &vm) final
Write from m_exp to output file.
Definition: OutputVtk.cpp:1613
void v_OutputFromPts(po::variables_map &vm) final
Write from pts to output file.
Definition: OutputVtk.cpp:1608
vtkUnstructuredGrid * GetVtkGrid()
Definition: OutputVtk.h:64
bool m_extraPlane
Flag if extra plane in case of fourier expansion in homogeneous dir.
Definition: OutputVtk.h:91
int m_numPlanes
Number of planes if homogeneous.
Definition: OutputVtk.h:88
std::string v_GetModuleName() final
Definition: OutputVtk.h:70
void v_OutputFromData(po::variables_map &vm) final
Write from data to output file.
Definition: OutputVtk.cpp:1601
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:1467
vtkSmartPointer< vtkUnstructuredGrid > OutputFromExpHighOrder(po::variables_map &vm)
Prepare high order Lagrange VTK output.
Definition: OutputVtk.cpp:1304
static ModuleKey m_className
Definition: OutputVtk.h:58
bool m_cachedMesh
Flag if mesh has been cached.
Definition: OutputVtk.h:94
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:1529
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:1026
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180