Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/dc/dca/VTKFile_8h_source.html
DOLFINx 0.6.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
VTKFile.h
1// Copyright (C) 2020 Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include <dolfinx/common/MPI.h>
10#include <filesystem>
11#include <functional>
12#include <memory>
13#include <string>
14
15namespace pugi
16{
17class xml_document;
18}
19
20namespace dolfinx::fem
21{
22template <typename T>
23class Function;
24}
25
26namespace dolfinx::mesh
27{
28class Mesh;
29}
30
31namespace dolfinx::io
32{
33
43{
44public:
46 VTKFile(MPI_Comm comm, const std::filesystem::path& filename,
47 const std::string& file_mode);
48
50 ~VTKFile();
51
53 void close();
54
56 void flush();
57
62 void write(const mesh::Mesh& mesh, double time = 0.0);
63
73 template <typename T>
74 void
75 write(const std::vector<std::reference_wrapper<const fem::Function<T>>>& u,
76 double t)
77 {
78 write_functions(u, t);
79 }
80
81private:
82 void write_functions(
83 const std::vector<std::reference_wrapper<const fem::Function<double>>>& u,
84 double t);
85 void write_functions(
86 const std::vector<
87 std::reference_wrapper<const fem::Function<std::complex<double>>>>& u,
88 double t);
89
90 std::unique_ptr<pugi::xml_document> _pvd_xml;
91
92 std::filesystem::path _filename;
93
94 // MPI communicator
96};
97} // namespace dolfinx::io
A duplicate MPI communicator and manage lifetime of the communicator.
Definition: MPI.h:42
This class represents a function in a finite element function space , given by.
Definition: Function.h:43
Output of meshes and functions in VTK/ParaView format. Isoparametric meshes of arbitrary degree are s...
Definition: VTKFile.h:43
void close()
Close file.
Definition: VTKFile.cpp:707
void write(const mesh::Mesh &mesh, double time=0.0)
Write a mesh to file. Supports arbitrary order Lagrange isoparametric cells.
Definition: VTKFile.cpp:737
~VTKFile()
Destructor.
Definition: VTKFile.cpp:697
void flush()
Flushes XML files to disk.
Definition: VTKFile.cpp:724
void write(const std::vector< std::reference_wrapper< const fem::Function< T > > > &u, double t)
Write finite elements function with an associated timestep.
Definition: VTKFile.h:75
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:34
Finite element method functionality.
Definition: assemble_matrix_impl.h:25
Support for file IO.
Definition: ADIOS2Writers.h:39
Mesh data structures and algorithms on meshes.
Definition: DofMap.h:31