Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
Nektar::FieldUtils::OutputVtk Class Referencefinal

Converter from fld to vtk. More...

#include <OutputVtk.h>

Inheritance diagram for Nektar::FieldUtils::OutputVtk:
[legend]

Public Member Functions

 OutputVtk (FieldSharedPtr f)
 
 ~OutputVtk () final=default
 
vtkUnstructuredGrid * GetVtkGrid ()
 
- Public Member Functions inherited from Nektar::FieldUtils::OutputVtkBase
 OutputVtkBase (FieldSharedPtr f)
 
 ~OutputVtkBase () override
 
- Public Member Functions inherited from Nektar::FieldUtils::OutputFileBase
 OutputFileBase (FieldSharedPtr f)
 
 ~OutputFileBase () override
 
- Public Member Functions inherited from Nektar::FieldUtils::OutputModule
 OutputModule (FieldSharedPtr p_f)
 
FIELD_UTILS_EXPORT void OpenStream ()
 Open a file for output. More...
 
- Public Member Functions inherited from Nektar::FieldUtils::Module
FIELD_UTILS_EXPORT Module (FieldSharedPtr p_f)
 
virtual ~Module ()=default
 
void Process (po::variables_map &vm)
 
std::string GetModuleName ()
 
std::string GetModuleDescription ()
 
const ConfigOptionGetConfigOption (const std::string &key) const
 
ModulePriority GetModulePriority ()
 
std::vector< ModuleKeyGetModulePrerequisites ()
 
FIELD_UTILS_EXPORT void RegisterConfig (std::string key, std::string value="")
 Register a configuration option with a module. More...
 
FIELD_UTILS_EXPORT void PrintConfig ()
 Print out all configuration options for a module. More...
 
FIELD_UTILS_EXPORT void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
FIELD_UTILS_EXPORT void AddFile (std::string fileType, std::string fileName)
 
FIELD_UTILS_EXPORT void EvaluateTriFieldAtEquiSpacedPts (LocalRegions::ExpansionSharedPtr &exp, const Array< OneD, const NekDouble > &infield, Array< OneD, NekDouble > &outfield)
 

Static Public Member Functions

static std::shared_ptr< Modulecreate (const FieldSharedPtr &f)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::FieldUtils::OutputVtkBase
static std::shared_ptr< Modulecreate (FieldSharedPtr f)
 Creates an instance of this class. More...
 

Static Public Attributes

static ModuleKey m_className
 
- Static Public Attributes inherited from Nektar::FieldUtils::OutputVtkBase
static ModuleKey m_className
 

Protected Member Functions

std::string v_GetModuleName () final
 
void v_OutputFromPts (po::variables_map &vm) final
 Write from pts to output file. More...
 
void v_OutputFromExp (po::variables_map &vm) final
 Write from m_exp to output file. More...
 
void v_OutputFromData (po::variables_map &vm) final
 Write from data to output file. More...
 
- Protected Member Functions inherited from Nektar::FieldUtils::OutputVtkBase
std::string v_GetModuleName () override
 
void v_OutputFromPts (po::variables_map &vm) override
 Write from pts to output file. More...
 
void v_OutputFromExp (po::variables_map &vm) override
 Write from m_exp to output file. More...
 
void v_OutputFromData (po::variables_map &vm) override
 Write from data to output file. More...
 
fs::path v_GetPath (std::string &filename, po::variables_map &vm) override
 
fs::path v_GetFullOutName (std::string &filename, po::variables_map &vm) override
 
std::string PrepareOutput (po::variables_map &vm)
 
- Protected Member Functions inherited from Nektar::FieldUtils::OutputFileBase
void v_Process (po::variables_map &vm) override
 Write fld to output file. More...
 
std::string v_GetModuleName () override
 
std::string v_GetModuleDescription () override
 
ModulePriority v_GetModulePriority () override
 
virtual void v_OutputFromPts (po::variables_map &vm)=0
 Write from pts to output file. More...
 
virtual void v_OutputFromExp (po::variables_map &vm)=0
 Write from m_exp to output file. More...
 
virtual void v_OutputFromData (po::variables_map &vm)=0
 Write from data to output file. More...
 
virtual fs::path v_GetPath (std::string &filename, po::variables_map &vm)
 
fs::path GetPath (std::string &filename, po::variables_map &vm)
 
virtual fs::path v_GetFullOutName (std::string &filename, po::variables_map &vm)
 
fs::path GetFullOutName (std::string &filename, po::variables_map &vm)
 
- Protected Member Functions inherited from Nektar::FieldUtils::Module
 Module ()
 
virtual void v_Process (po::variables_map &vm)
 
virtual std::string v_GetModuleName ()
 
virtual std::string v_GetModuleDescription ()
 
virtual ModulePriority v_GetModulePriority ()
 
virtual std::vector< ModuleKeyv_GetModulePrerequisites ()
 

Protected Attributes

vtkSmartPointer< vtkUnstructuredGrid > m_vtkMesh
 Cache file for unstructured grid VTK mesh data. More...
 
int m_numPlanes = 1
 Number of planes if homogeneous. More...
 
bool m_extraPlane = false
 Flag if extra plane in case of fourier expansion in homogeneous dir. More...
 
bool m_cachedMesh = false
 Flag if mesh has been cached. More...
 
- Protected Attributes inherited from Nektar::FieldUtils::OutputFileBase
bool m_requireEquiSpaced
 
bool m_prohibitWrite = false
 
bool m_equispacedSetup = false
 
- Protected Attributes inherited from Nektar::FieldUtils::OutputModule
std::ofstream m_fldFile
 Output stream. More...
 
- Protected Attributes inherited from Nektar::FieldUtils::Module
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 
std::set< std::string > m_allowedFiles
 List of allowed file formats. More...
 

Private Member Functions

vtkSmartPointer< vtkUnstructuredGrid > OutputFromExpHighOrder (po::variables_map &vm)
 Prepare high order Lagrange VTK output. More...
 
void AddFieldDataToVTKHighOrder (po::variables_map &vm, std::string &filename, vtkSmartPointer< vtkUnstructuredGrid > &vtkMesh)
 Add field data to high order Lagrange VTK output. More...
 
vtkSmartPointer< vtkUnstructuredGrid > OutputFromExpLowOrder ()
 Prepare low order VTK output. More...
 
void AddFieldDataToVTKLowOrder (po::variables_map &vm, std::string &filename, vtkSmartPointer< vtkUnstructuredGrid > &vtkMesh)
 Add field data to low order VTK output. More...
 
void OutputFromExpLowOrderMultiBlock (po::variables_map &vm, std::string &filename)
 Prepare low order multi-block VTK output & add field data. More...
 
void WriteVTK (vtkDataObject *vtkMesh, std::string &filename, po::variables_map &vm)
 Write VTK file using. More...
 
void WritePVtu (po::variables_map &vm)
 Write the parallel .pvtu file. More...
 

Additional Inherited Members

- Public Attributes inherited from Nektar::FieldUtils::Module
FieldSharedPtr m_f
 Field object. More...
 

Detailed Description

Converter from fld to vtk.

Definition at line 49 of file OutputVtk.h.

Constructor & Destructor Documentation

◆ OutputVtk()

Nektar::FieldUtils::OutputVtk::OutputVtk ( FieldSharedPtr  f)
explicit

Definition at line 56 of file OutputVtk.cpp.

56 : OutputVtkBase(std::move(f))
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}
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272

References Nektar::FieldUtils::Module::m_config, and Nektar::FieldUtils::OutputFileBase::m_requireEquiSpaced.

◆ ~OutputVtk()

Nektar::FieldUtils::OutputVtk::~OutputVtk ( )
finaldefault

Member Function Documentation

◆ AddFieldDataToVTKHighOrder()

void Nektar::FieldUtils::OutputVtk::AddFieldDataToVTKHighOrder ( po::variables_map &  vm,
std::string &  filename,
vtkSmartPointer< vtkUnstructuredGrid > &  vtkMesh 
)
private

Add field data to high order Lagrange VTK output.

Definition at line 1433 of file OutputVtk.cpp.

1436{
1437 // Convert expansion to an equispaced points field
1438 if (m_cachedMesh)
1439 {
1440 auto equispaced =
1442 equispaced->Process(vm);
1443 }
1444 LibUtilities::PtsFieldSharedPtr fPts = m_f->m_fieldPts;
1445
1446 int nPts = static_cast<int>(fPts->GetNpoints());
1447 int dim = static_cast<int>(fPts->GetDim());
1448
1449 Array<OneD, Array<OneD, NekDouble>> pts;
1450 fPts->GetPts(pts);
1451
1452 // Insert field information
1453 for (int i = 0; i < fPts->GetNFields(); ++i)
1454 {
1455 vtkNew<vtkDoubleArray> fieldData;
1456 fieldData->SetArray(&pts[dim + i][0], nPts, 1);
1457 fieldData->SetName(&fPts->GetFieldName(i)[0]);
1458 vtkMesh->GetPointData()->AddArray(fieldData.GetPointer());
1459 }
1460
1461 if (!m_prohibitWrite)
1462 {
1463 WriteVTK(vtkMesh, filename, vm);
1464 }
1465}
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
void WriteVTK(vtkDataObject *vtkMesh, std::string &filename, po::variables_map &vm)
Write VTK file using.
Definition: OutputVtk.cpp:1467
bool m_cachedMesh
Flag if mesh has been cached.
Definition: OutputVtk.h:94
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:184

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), m_cachedMesh, Nektar::FieldUtils::Module::m_f, Nektar::FieldUtils::OutputFileBase::m_prohibitWrite, and WriteVTK().

Referenced by v_OutputFromExp().

◆ AddFieldDataToVTKLowOrder()

void Nektar::FieldUtils::OutputVtk::AddFieldDataToVTKLowOrder ( po::variables_map &  vm,
std::string &  filename,
vtkSmartPointer< vtkUnstructuredGrid > &  vtkMesh 
)
private

Add field data to low order VTK output.

Definition at line 975 of file OutputVtk.cpp.

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 if (!m_prohibitWrite)
1017 {
1018 WriteVTK(vtkMesh, filename, vm);
1019 }
1020}
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

References m_extraPlane, Nektar::FieldUtils::Module::m_f, m_numPlanes, Nektar::FieldUtils::OutputFileBase::m_prohibitWrite, and WriteVTK().

Referenced by v_OutputFromExp().

◆ create()

static std::shared_ptr< Module > Nektar::FieldUtils::OutputVtk::create ( const FieldSharedPtr f)
inlinestatic

Creates an instance of this class.

Definition at line 53 of file OutputVtk.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ GetVtkGrid()

vtkUnstructuredGrid * Nektar::FieldUtils::OutputVtk::GetVtkGrid ( )
inline

Definition at line 64 of file OutputVtk.h.

65 {
66 return m_vtkMesh.GetPointer();
67 }
vtkSmartPointer< vtkUnstructuredGrid > m_vtkMesh
Cache file for unstructured grid VTK mesh data.
Definition: OutputVtk.h:85

References m_vtkMesh.

◆ OutputFromExpHighOrder()

vtkSmartPointer< vtkUnstructuredGrid > Nektar::FieldUtils::OutputVtk::OutputFromExpHighOrder ( po::variables_map &  vm)
private

Prepare high order Lagrange VTK output.

Definition at line 1304 of file OutputVtk.cpp.

1306{
1307 ASSERTL0(
1308 m_f->m_numHomogeneousDir == 0,
1309 "High order VTK is not implemented for homogeneous expansion types.")
1310
1311 // Save shapetypes before convert to equispaced because that nukes the
1312 // explist
1313 std::vector<LibUtilities::ShapeType> sType;
1314 for (auto &i : *m_f->m_exp[0]->GetExp())
1315 {
1316 sType.emplace_back(i->GetGeom()->GetShapeType());
1317 }
1318
1319 // Convert expansion to an equispaced points field
1320 auto equispaced =
1322 equispaced->Process(vm);
1323
1324 // Insert points
1325 vtkNew<vtkPoints> vtkPoints;
1326 LibUtilities::PtsFieldSharedPtr fPts = m_f->m_fieldPts;
1327 vtkPoints->SetDataType(VTK_DOUBLE);
1328
1329 Array<OneD, Array<OneD, NekDouble>> pts;
1330 fPts->GetPts(pts);
1331
1332 int dim = static_cast<int>(fPts->GetDim());
1333 int nPts = static_cast<int>(fPts->GetNpoints());
1334 switch (dim)
1335 {
1336 case 3:
1337 for (int i = 0; i < nPts; ++i)
1338 {
1339 vtkPoints->InsertNextPoint(pts[0][i], pts[1][i], pts[2][i]);
1340 }
1341 break;
1342 case 2:
1343 for (int i = 0; i < nPts; ++i)
1344 {
1345 vtkPoints->InsertNextPoint(pts[0][i], pts[1][i], 0.0);
1346 }
1347 break;
1348 case 1:
1349 for (int i = 0; i < nPts; ++i)
1350 {
1351 vtkPoints->InsertNextPoint(pts[0][i], 0.0, 0.0);
1352 }
1353 break;
1354 default:
1355 break;
1356 }
1357
1358 // Mesh file to output
1359 vtkSmartPointer<vtkUnstructuredGrid> vtkMesh =
1360 vtkSmartPointer<vtkUnstructuredGrid>::New();
1361 vtkMesh->SetPoints(vtkPoints.GetPointer());
1362
1363 // Cache ordering for shape type & npts so we aren't recreating mappings
1364 std::unordered_map<size_t, std::vector<long long>> mappingCache;
1365
1366 // Get offset per element in a vector from ppe
1367 std::vector<int> ppe = m_f->m_fieldPts->GetPointsPerElement();
1368 Array<OneD, int> ppeOffset(ppe.size() + 1, 0.0);
1369 std::partial_sum(ppe.begin(), ppe.end(), &ppeOffset[1]);
1370
1371 // Insert elements
1372 for (int i = 0; i < ppe.size(); ++i)
1373 {
1374 // Construct inverse of input reordering.
1375 // First try to find it in our mapping cache.
1376 auto oIt = mappingCache.find(key2(sType[i], ppe[i]));
1377 if (oIt == mappingCache.end())
1378 {
1379 std::vector<long long> p(ppe[i]);
1380 std::iota(p.begin(), p.end(), 0);
1381 switch (sType[i])
1382 {
1384 p = quadTensorNodeOrdering(p);
1385 break;
1387 p = triTensorNodeOrdering(p);
1388 break;
1390 p = tetTensorNodeOrdering(p);
1391 break;
1393 p = hexTensorNodeOrdering(p);
1394 break;
1396 p = prismTensorNodeOrdering(p);
1397 break;
1398 default:
1400 "High-order VTU output not set up for the " +
1401 static_cast<std::string>(
1402 LibUtilities::ShapeTypeMap[sType[i]]) +
1403 " shape type.");
1404 break;
1405 }
1406
1407 // Invert the ordering as this is for spectral -> recursive (VTU)
1408 std::vector<long long> inv(ppe[i]);
1409 for (int j = 0; j < ppe[i]; ++j)
1410 {
1411 inv[p[j]] = j;
1412 }
1413
1414 oIt =
1415 mappingCache.insert(std::make_pair(key2(sType[i], ppe[i]), inv))
1416 .first;
1417 }
1418
1419 // Add offset to reordering
1420 std::vector<long long> p = oIt->second;
1421 for (long long &j : p)
1422 {
1423 j += ppeOffset[i];
1424 }
1425
1426 vtkMesh->InsertNextCell(GetHighOrderVtkCellType(sType[i]), ppe[i],
1427 &p[0]);
1428 }
1429
1430 return vtkMesh;
1431}
#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
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:75
STL namespace.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, Nektar::FieldUtils::Module::m_f, NEKERROR, CellMLToNektar.cellml_metadata::p, and Nektar::LibUtilities::ShapeTypeMap.

Referenced by v_OutputFromExp().

◆ OutputFromExpLowOrder()

vtkSmartPointer< vtkUnstructuredGrid > Nektar::FieldUtils::OutputVtk::OutputFromExpLowOrder ( )
private

Prepare low order VTK output.

Definition at line 852 of file OutputVtk.cpp.

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() ==
882 LibUtilities::eFourier &&
883 m_f->m_exp[0]->GetHomogeneousBasis()->GetPointsType() ==
884 LibUtilities::eFourierEvenlySpaced);
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
929 Array<OneD, NekDouble> tmp;
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}
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
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

References ASSERTL0, Nektar::UnitTests::d(), Nektar::MultiRegions::e3DH1D, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierEvenlySpaced, Vmath::Fill(), m_extraPlane, Nektar::FieldUtils::Module::m_f, m_numPlanes, CellMLToNektar.cellml_metadata::p, Vmath::Vcopy(), and Nektar::UnitTests::z().

Referenced by v_OutputFromExp().

◆ OutputFromExpLowOrderMultiBlock()

void Nektar::FieldUtils::OutputVtk::OutputFromExpLowOrderMultiBlock ( po::variables_map &  vm,
std::string &  filename 
)
private

Prepare low order multi-block VTK output & add field data.

Definition at line 1022 of file OutputVtk.cpp.

1024{
1025 ASSERTL0(
1026 m_f->m_numHomogeneousDir == 0,
1027 "Multi block VTK is not implemented for homogeneous expansion types.")
1028
1029 ASSERTL0(m_f->m_comm->GetSpaceComm()->IsSerial(),
1030 "Multi block VTK is not implemented in parallel.")
1031
1032 int dim = m_f->m_graph->GetMeshDimension();
1033
1034 // Create mappings from geometry id to expansion ids
1035 std::array<std::map<int, std::pair<int, int>>, 4> geomIdToExpId;
1036 int nElmts = static_cast<int>(m_f->m_exp[0]->GetNumElmts());
1037 for (int i = 0; i < nElmts; ++i)
1038 {
1039 auto geom = m_f->m_exp[0]->GetExp(i)->GetGeom();
1040 geomIdToExpId[geom->GetShapeDim()][geom->GetGlobalID()] =
1041 std::make_pair(i, -1);
1042
1043 for (int j = 0; j < geom->GetNumFaces(); ++j)
1044 {
1045 geomIdToExpId[2][geom->GetFid(j)] = std::make_pair(i, j);
1046 }
1047
1048 for (int j = 0; j < geom->GetNumEdges(); ++j)
1049 {
1050 geomIdToExpId[1][geom->GetEid(j)] = std::make_pair(i, j);
1051 }
1052 }
1053
1054 // Cache ordering so we aren't recreating mappings
1055 std::unordered_map<size_t, std::vector<long long>> mappingCache;
1056
1057 std::map<int, vtkNew<vtkUnstructuredGrid>> vtkMesh;
1058 std::map<int, SpatialDomains::CompositeSharedPtr> composites =
1059 m_f->m_graph->GetComposites();
1060 std::map<int, std::string> compositeNames;
1061 for (auto &comp : composites)
1062 {
1063 // Vector of field data
1064 std::vector<vtkNew<vtkDoubleArray>> fieldData(m_f->m_variables.size());
1065
1066 // Insert points
1067 vtkNew<vtkPoints> vtkPoints;
1068 vtkPoints->SetDataType(VTK_DOUBLE);
1069
1070 int compId = comp.first;
1071 std::vector<std::shared_ptr<SpatialDomains::Geometry>> geomVec =
1072 comp.second->m_geomVec;
1073
1074 unsigned int offset = 0;
1075 for (auto &geom : geomVec)
1076 {
1077 int geomId = geom->GetGlobalID();
1078 int sDim = geom->GetShapeDim();
1079 auto type = (sDim == 3) ? VTK_HEXAHEDRON
1080 : (sDim == 2) ? VTK_QUAD
1081 : VTK_LINE;
1082 int nQpts = (sDim == 3) ? 8 : (sDim == 2) ? 4 : 2;
1083
1085 m_f->m_exp[0]->GetExp(geomIdToExpId[sDim][geomId].first);
1086
1087 unsigned int nPts = exp->GetTotPoints();
1088 Array<OneD, NekDouble> x(nPts, 0.0), y(nPts, 0.0), z(nPts, 0.0);
1089 exp->GetCoords(x, y, z);
1090
1091 int offsetPhys = m_f->m_exp[0]->GetPhys_Offset(
1092 geomIdToExpId[sDim][geomId].first);
1093 if (sDim == dim)
1094 {
1095 for (int j = 0; j < nPts; ++j)
1096 {
1097 vtkPoints->InsertNextPoint(x[j], y[j], z[j]);
1098
1099 // Add field data
1100 for (int k = 0; k < m_f->m_variables.size(); ++k)
1101 {
1102 fieldData[k]->InsertNextValue(
1103 m_f->m_exp[k]->GetPhys()[j + offsetPhys]);
1104 }
1105 }
1106 }
1107 else
1108 {
1109 Array<OneD, int> pointsMap;
1110 exp->GetTracePhysMap(geomIdToExpId[sDim][geomId].second,
1111 pointsMap);
1112 for (int j : pointsMap)
1113 {
1114 vtkPoints->InsertNextPoint(x[j], y[j], z[j]);
1115
1116 // Add field data
1117 for (int k = 0; k < m_f->m_variables.size(); ++k)
1118 {
1119 fieldData[k]->InsertNextValue(
1120 m_f->m_exp[k]->GetPhys()[offsetPhys + j]);
1121 }
1122 }
1123
1124 exp = exp->GetTraceExp(geomIdToExpId[sDim][geomId].second);
1125 nPts = pointsMap.size();
1126 }
1127
1128 Array<OneD, int> nquad(3, 0);
1129 for (int j = 0; j < sDim; ++j)
1130 {
1131 nquad[j] = exp->GetNumPoints(j);
1132 }
1133
1134 auto oIt = mappingCache.find(key3(nquad[0], nquad[1], nquad[2]));
1135 if (oIt == mappingCache.end())
1136 {
1137 auto p = lowOrderMapping(sDim, nquad);
1138 oIt = mappingCache
1139 .insert(std::make_pair(
1140 key3(nquad[0], nquad[1], nquad[2]), p))
1141 .first;
1142 }
1143
1144 auto p = oIt->second;
1145 std::for_each(p.begin(), p.end(),
1146 [offset](long long &d) { d += offset; });
1147
1148 for (int j = 0; j < p.size(); j += nQpts)
1149 {
1150 vtkMesh[compId]->InsertNextCell(type, nQpts, &p[j]);
1151 }
1152
1153 offset += nPts;
1154 }
1155
1156 vtkMesh[compId]->SetPoints(vtkPoints.GetPointer());
1157
1158 // Add all fields to vtkMesh
1159 for (int k = 0; k < m_f->m_variables.size(); ++k)
1160 {
1161 fieldData[k]->SetName(&m_f->m_variables[k][0]);
1162 vtkMesh[compId]->GetPointData()->AddArray(
1163 fieldData[k].GetPointer());
1164 }
1165
1166 // Find composite names if exists and store in map
1167 // Set name as boundary label if it exists otherwise index ID
1168 compositeNames[compId] = "Composite ID " + std::to_string(compId);
1169 auto clabels = m_f->m_graph->GetCompositesLabels();
1170 auto oIt = clabels.find(compId);
1171 if (oIt != clabels.end())
1172 {
1173 compositeNames[compId] = oIt->second;
1174 }
1175 }
1176
1177 // Main multi-block set
1178 vtkNew<vtkMultiBlockDataSet> mainBlock;
1179
1180 std::set<int> compSet;
1181
1182 // Fill domain multi blocks from composites
1183 vtkNew<vtkMultiBlockDataSet> mainDomainBlock;
1184 auto domains = m_f->m_graph->GetDomain();
1185 std::vector<vtkNew<vtkMultiBlockDataSet>> domainMultiBlocks(domains.size());
1186 for (int i = 0; i < domains.size(); ++i)
1187 {
1188 auto dom = domains[i];
1189
1190 // Loop over composites and see if in domain
1191 for (auto &comp : composites)
1192 {
1193 int compId = comp.first;
1194 if (dom.find(compId) != dom.end())
1195 {
1196 unsigned int nBlock = domainMultiBlocks[i]->GetNumberOfBlocks();
1197 domainMultiBlocks[i]->SetBlock(nBlock,
1198 vtkMesh[compId].GetPointer());
1199 domainMultiBlocks[i]->GetMetaData(nBlock)->Set(
1200 vtkCompositeDataSet::NAME(),
1201 compositeNames[compId].c_str());
1202 compSet.insert(compId);
1203 }
1204 }
1205
1206 unsigned int nBlock = mainDomainBlock->GetNumberOfBlocks();
1207 mainDomainBlock->SetBlock(nBlock, domainMultiBlocks[i].GetPointer());
1208 mainDomainBlock->GetMetaData(nBlock)->Set(
1209 vtkCompositeDataSet::NAME(),
1210 ("Domain ID " + std::to_string(i)).c_str());
1211 }
1212
1213 if (mainDomainBlock->GetNumberOfBlocks() != 0)
1214 {
1215 auto nBlock = static_cast<unsigned int>(mainBlock->GetNumberOfBlocks());
1216 mainBlock->SetBlock(nBlock, mainDomainBlock.GetPointer());
1217 mainBlock->GetMetaData(nBlock)->Set(vtkCompositeDataSet::NAME(),
1218 "Domains");
1219 }
1220
1221 // Fill boundary multi blocks from composites
1222 SpatialDomains::BoundaryConditions bcs(m_f->m_session,
1223 m_f->m_exp[0]->GetGraph());
1225 bcs.GetBoundaryRegions();
1226
1227 vtkNew<vtkMultiBlockDataSet> mainBoundaryBlock;
1228 std::vector<vtkNew<vtkMultiBlockDataSet>> boundaryMultiBlocks(
1229 bregions.size());
1230 int cnt = 0;
1231 for (auto &boundary : bregions)
1232 {
1233 // Loop over composites and see if in boundary
1234 for (auto &comp : composites)
1235 {
1236 int compId = comp.first;
1237 if (boundary.second->find(compId) != boundary.second->end())
1238 {
1239 unsigned int nBlock =
1240 boundaryMultiBlocks[cnt]->GetNumberOfBlocks();
1241 boundaryMultiBlocks[cnt]->SetBlock(
1242 nBlock, vtkMesh[compId].GetPointer());
1243 boundaryMultiBlocks[cnt]->GetMetaData(nBlock)->Set(
1244 vtkCompositeDataSet::NAME(),
1245 compositeNames[compId].c_str());
1246 compSet.insert(compId);
1247 }
1248 }
1249
1250 unsigned int nBlock = mainBoundaryBlock->GetNumberOfBlocks();
1251 mainBoundaryBlock->SetBlock(nBlock,
1252 boundaryMultiBlocks[cnt++].GetPointer());
1253
1254 // Set name as boundary label if it exists otherwise index ID
1255 std::string name = "Boundary ID " + std::to_string(boundary.first);
1256 auto blabels = bcs.GetBoundaryLabels();
1257 auto oIt = blabels.find(boundary.first);
1258 if (oIt != blabels.end())
1259 {
1260 name = oIt->second;
1261 }
1262
1263 mainBoundaryBlock->GetMetaData(nBlock)->Set(vtkCompositeDataSet::NAME(),
1264 name.c_str());
1265 }
1266
1267 if (mainBoundaryBlock->GetNumberOfBlocks() != 0)
1268 {
1269 auto nBlock = static_cast<unsigned int>(mainBlock->GetNumberOfBlocks());
1270 mainBlock->SetBlock(nBlock, mainBoundaryBlock.GetPointer());
1271 mainBlock->GetMetaData(nBlock)->Set(vtkCompositeDataSet::NAME(),
1272 "Boundaries");
1273 }
1274
1275 // Include all other composites (not domains or boundaries)
1276 vtkNew<vtkMultiBlockDataSet> compsMultiBlocks;
1277 cnt = 0;
1278 for (auto &comp : composites)
1279 {
1280 int compId = comp.first;
1281 if (compSet.find(compId) == compSet.end())
1282 {
1283 unsigned int nBlock = compsMultiBlocks->GetNumberOfBlocks();
1284 compsMultiBlocks->SetBlock(nBlock, vtkMesh[compId].GetPointer());
1285 compsMultiBlocks->GetMetaData(nBlock)->Set(
1286 vtkCompositeDataSet::NAME(), compositeNames[compId].c_str());
1287 }
1288 }
1289
1290 if (compsMultiBlocks->GetNumberOfBlocks() != 0)
1291 {
1292 auto nBlock = static_cast<unsigned int>(mainBlock->GetNumberOfBlocks());
1293 mainBlock->SetBlock(nBlock, compsMultiBlocks.GetPointer());
1294 mainBlock->GetMetaData(nBlock)->Set(vtkCompositeDataSet::NAME(),
1295 "Other composites");
1296 }
1297
1298 if (!m_prohibitWrite)
1299 {
1300 WriteVTK(mainBlock.GetPointer(), filename, vm);
1301 }
1302}
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:210

References ASSERTL0, Nektar::UnitTests::d(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryLabels(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), Nektar::FieldUtils::Module::m_f, Nektar::FieldUtils::OutputFileBase::m_prohibitWrite, CellMLToNektar.pycml::name, CellMLToNektar.cellml_metadata::p, WriteVTK(), and Nektar::UnitTests::z().

Referenced by v_OutputFromExp().

◆ v_GetModuleName()

std::string Nektar::FieldUtils::OutputVtk::v_GetModuleName ( )
inlinefinalprotectedvirtual

Reimplemented from Nektar::FieldUtils::OutputFileBase.

Definition at line 70 of file OutputVtk.h.

71 {
72 return "OutputVtk";
73 }

◆ v_OutputFromData()

void Nektar::FieldUtils::OutputVtk::v_OutputFromData ( po::variables_map &  vm)
finalprotectedvirtual

Write from data to output file.

Implements Nektar::FieldUtils::OutputFileBase.

Definition at line 1601 of file OutputVtk.cpp.

1602{
1604 "OutputVtk can't write using only FieldData. You may need "
1605 "to add a mesh XML file to your input files.");
1606}

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_OutputFromExp()

void Nektar::FieldUtils::OutputVtk::v_OutputFromExp ( po::variables_map &  vm)
finalprotectedvirtual

Write from m_exp to output file.

Implements Nektar::FieldUtils::OutputFileBase.

Definition at line 1613 of file OutputVtk.cpp.

1614{
1615 // Move geometry based on zones data in .xml and time in .fld metadatamap
1616 // Perform movement of zones based on time in field files metadata map
1617 if (m_f->m_graph->GetMovement()->GetMoveFlag())
1618 {
1619 if (!m_f->m_fieldMetaDataMap["Time"].empty())
1620 {
1621 m_f->m_graph->GetMovement()->PerformMovement(
1622 std::stod(m_f->m_fieldMetaDataMap["Time"]));
1623 for (auto &i : m_f->m_exp)
1624 {
1625 i->Reset();
1626 }
1627 }
1628 }
1629
1630 if (m_config["legacy"].m_beenSet)
1631 {
1632 ASSERTL0(!m_config["multiblock"].m_beenSet,
1633 "Multi block VTK is not implemented for legacy output.")
1634
1635 ASSERTL0(!m_config["highorder"].m_beenSet,
1636 "High order VTK is not implemented for legacy output.")
1637
1638 // No caching of mesh data in legacy output
1640 return;
1641 }
1642
1643 std::string filename;
1644
1645 // Extract the output filename and extension
1646 if (!m_prohibitWrite)
1647 {
1648 filename = OutputVtkBase::PrepareOutput(vm);
1649 }
1650
1651 // Save mesh state (if using filter this allows us to only ProcessEquispaced
1652 // if needed)
1653 if (!m_f->m_graph->GetMovement()->GetMoveFlag())
1654 {
1656 }
1657
1658 if (m_config["highorder"].m_beenSet)
1659 {
1660 ASSERTL0(!m_config["multiblock"].m_beenSet,
1661 "Multi block VTK is not implemented for high-order output.")
1662
1663 if (!m_cachedMesh)
1664 {
1666 }
1667
1669 }
1670 else if (m_config["multiblock"].m_beenSet)
1671 {
1672 // No mesh caching for multiblock due to complexity of field data
1673 OutputFromExpLowOrderMultiBlock(vm, filename);
1674 }
1675 else
1676 {
1677 if (!m_cachedMesh)
1678 {
1680 }
1681
1682 AddFieldDataToVTKLowOrder(vm, filename, m_vtkMesh);
1683 }
1684}
std::string PrepareOutput(po::variables_map &vm)
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
void v_OutputFromExp(po::variables_map &vm) final
Write from m_exp to output file.
Definition: OutputVtk.cpp:1613
void AddFieldDataToVTKLowOrder(po::variables_map &vm, std::string &filename, vtkSmartPointer< vtkUnstructuredGrid > &vtkMesh)
Add field data to low order VTK output.
Definition: OutputVtk.cpp:975
vtkSmartPointer< vtkUnstructuredGrid > OutputFromExpHighOrder(po::variables_map &vm)
Prepare high order Lagrange VTK output.
Definition: OutputVtk.cpp:1304

References AddFieldDataToVTKHighOrder(), AddFieldDataToVTKLowOrder(), ASSERTL0, m_cachedMesh, Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, Nektar::FieldUtils::OutputFileBase::m_prohibitWrite, m_vtkMesh, OutputFromExpHighOrder(), OutputFromExpLowOrder(), OutputFromExpLowOrderMultiBlock(), Nektar::FieldUtils::OutputVtkBase::PrepareOutput(), and Nektar::FieldUtils::OutputVtkBase::v_OutputFromExp().

◆ v_OutputFromPts()

void Nektar::FieldUtils::OutputVtk::v_OutputFromPts ( po::variables_map &  vm)
finalprotectedvirtual

Write from pts to output file.

Implements Nektar::FieldUtils::OutputFileBase.

Definition at line 1608 of file OutputVtk.cpp.

1609{
1611}
void v_OutputFromPts(po::variables_map &vm) override
Write from pts to output file.

References Nektar::FieldUtils::OutputVtkBase::v_OutputFromPts().

◆ WritePVtu()

void Nektar::FieldUtils::OutputVtk::WritePVtu ( po::variables_map &  vm)
private

Write the parallel .pvtu file.

Definition at line 1529 of file OutputVtk.cpp.

1530{
1531 std::string filename = m_config["outfile"].as<std::string>();
1532 int dot = static_cast<int>(filename.find_last_of('.'));
1533 std::string body = filename.substr(0, dot);
1534 filename = body + ".pvtu";
1535
1536 std::ofstream outfile(filename.c_str());
1537
1538 int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
1539 std::string path =
1541
1542 outfile << "<?xml version=\"1.0\"?>" << endl;
1543 outfile << R"(<VTKFile type="PUnstructuredGrid" version="0.1" )"
1544 << "byte_order=\"LittleEndian\">" << endl;
1545 outfile << " <PUnstructuredGrid GhostLevel=\"0\">" << endl;
1546
1547 // Add time if exists
1548 if (!m_f->m_fieldMetaDataMap["Time"].empty())
1549 {
1550 outfile << " <FieldData> " << endl;
1551 outfile << " <DataArray type=\"Float64\" Name=\"TimeValue\" "
1552 "NumberOfTuples=\"1\" format=\"ascii\">"
1553 << endl;
1554 outfile << m_f->m_fieldMetaDataMap["Time"] << std::endl;
1555 outfile << " </DataArray>" << std::endl;
1556 outfile << " </FieldData> " << endl;
1557 }
1558
1559 // Add point coordinates
1560 outfile << " <PPoints> " << endl;
1561 outfile << R"( <PDataArray type="Float64" NumberOfComponents=")" << 3
1562 << "\"/> " << endl;
1563 outfile << " </PPoints>" << endl;
1564
1565 // Add cell information
1566 outfile << " <PCells>" << endl;
1567 outfile << " <PDataArray type=\"Int64\" Name=\"connectivity\" "
1568 "NumberOfComponents=\"1\"/>"
1569 << endl;
1570 outfile << " <PDataArray type=\"Int64\" Name=\"offsets\" "
1571 "NumberOfComponents=\"1\"/>"
1572 << endl;
1573 outfile << " <PDataArray type=\"UInt8\" Name=\"types\" "
1574 "NumberOfComponents=\"1\"/>"
1575 << endl;
1576 outfile << " </PCells>" << endl;
1577
1578 // Add fields (point data)
1579 outfile << " <PPointData>" << endl;
1580 for (auto &var : m_f->m_variables)
1581 {
1582 outfile << R"( <PDataArray type="Float64" Name=")" << var << "\"/>"
1583 << endl;
1584 }
1585 outfile << " </PPointData>" << endl;
1586
1587 // Add parallel files
1588 for (int i = 0; i < nprocs; ++i)
1589 {
1590 boost::format pad("P%1$07d.vtu");
1591 pad % i;
1592 outfile << " <Piece Source=\"" << path << "/" << pad.str() << "\"/>"
1593 << endl;
1594 }
1595 outfile << " </PUnstructuredGrid>" << endl;
1596 outfile << "</VTKFile>" << endl;
1597
1598 cout << "Written file: " << filename << endl;
1599}
fs::path v_GetPath(std::string &filename, po::variables_map &vm) override
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
Definition: Filesystem.hpp:56

References CellMLToNektar.pycml::format, Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, Nektar::LibUtilities::PortablePath(), and Nektar::FieldUtils::OutputVtkBase::v_GetPath().

Referenced by WriteVTK().

◆ WriteVTK()

void Nektar::FieldUtils::OutputVtk::WriteVTK ( vtkDataObject *  vtkMesh,
std::string &  filename,
po::variables_map &  vm 
)
private

Write VTK file using.

Parameters
vtkMesh

Definition at line 1467 of file OutputVtk.cpp.

1469{
1470 // Initialise base writer class as default structured grid writer
1471 vtkSmartPointer<vtkXMLWriter> vtkMeshWriter =
1472 vtkNew<vtkXMLUnstructuredGridWriter>().GetPointer();
1473
1474 // If multi block we need to switch to the MultiBlock writer
1475 if (m_config["multiblock"].m_beenSet)
1476 {
1477 vtkMeshWriter = vtkNew<vtkXMLMultiBlockDataWriter>().GetPointer();
1478 }
1479
1480 // Choose the correct extension based on which mesh writer is used
1481 int dot = static_cast<int>(filename.find_last_of('.'));
1482 std::string body = filename.substr(0, dot + 1);
1483 filename = body + vtkMeshWriter->GetDefaultFileExtension();
1484
1485 // Set time information if exists
1486 if (!m_f->m_fieldMetaDataMap["Time"].empty())
1487 {
1488 vtkMesh->GetInformation()->Set(
1489 vtkDataObject::DATA_TIME_STEP(),
1490 std::stod(m_f->m_fieldMetaDataMap["Time"]));
1491 }
1492 else if (!m_f->m_fieldMetaDataMap["FinalTime"].empty())
1493 {
1494 vtkMesh->GetInformation()->Set(
1495 vtkDataObject::DATA_TIME_STEP(),
1496 std::stod(m_f->m_fieldMetaDataMap["FinalTime"]));
1497 }
1498
1499 // Write out the new mesh in XML format (don't support legacy
1500 // format here as we still have standard OutputVtk.cpp)
1501 vtkMeshWriter->SetFileName(filename.c_str());
1502
1503#if VTK_MAJOR_VERSION <= 5
1504 vtkMeshWriter->SetInput(vtkMesh);
1505#else
1506 vtkMeshWriter->SetInputData(vtkMesh);
1507#endif
1508
1509 if (m_config["uncompress"].m_beenSet)
1510 {
1511 vtkMeshWriter->SetDataModeToAscii();
1512 }
1513
1514 vtkMeshWriter->Update();
1515
1516 std::cout << "Written file: " << filename << std::endl;
1517
1518 // Output parallel file information if needed - using a manual write.
1519 // We could use the VTK lib to do this, but that requires VTK with MPI
1520 // enabled & messing about with the parallel controller & changing
1521 // our file naming scheme as VTK forces _${proc-number} as a suffix...
1522 if (m_f->m_comm->GetSpaceComm()->TreatAsRankZero() &&
1523 !m_f->m_comm->GetSpaceComm()->IsSerial())
1524 {
1525 WritePVtu(vm);
1526 }
1527}
void WritePVtu(po::variables_map &vm)
Write the parallel .pvtu file.
Definition: OutputVtk.cpp:1529

References Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, and WritePVtu().

Referenced by AddFieldDataToVTKHighOrder(), AddFieldDataToVTKLowOrder(), and OutputFromExpLowOrderMultiBlock().

Member Data Documentation

◆ m_cachedMesh

bool Nektar::FieldUtils::OutputVtk::m_cachedMesh = false
protected

Flag if mesh has been cached.

Definition at line 94 of file OutputVtk.h.

Referenced by AddFieldDataToVTKHighOrder(), and v_OutputFromExp().

◆ m_className

ModuleKey Nektar::FieldUtils::OutputVtk::m_className
static
Initial value:
ModuleKey(eOutputModule, "vtu"), OutputVtk::create, "Writes a VTU file.")
static std::shared_ptr< Module > create(const FieldSharedPtr &f)
Creates an instance of this class.
Definition: OutputVtk.h:53
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47

Definition at line 58 of file OutputVtk.h.

◆ m_extraPlane

bool Nektar::FieldUtils::OutputVtk::m_extraPlane = false
protected

Flag if extra plane in case of fourier expansion in homogeneous dir.

Definition at line 91 of file OutputVtk.h.

Referenced by AddFieldDataToVTKLowOrder(), and OutputFromExpLowOrder().

◆ m_numPlanes

int Nektar::FieldUtils::OutputVtk::m_numPlanes = 1
protected

Number of planes if homogeneous.

Definition at line 88 of file OutputVtk.h.

Referenced by AddFieldDataToVTKLowOrder(), and OutputFromExpLowOrder().

◆ m_vtkMesh

vtkSmartPointer<vtkUnstructuredGrid> Nektar::FieldUtils::OutputVtk::m_vtkMesh
protected

Cache file for unstructured grid VTK mesh data.

Definition at line 85 of file OutputVtk.h.

Referenced by GetVtkGrid(), and v_OutputFromExp().