OpenLB 1.8.1
Loading...
Searching...
No Matches
vtuSurfaceWriter.h
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2024 Christoph Gaul, Adrian Kummerlaender
4 * E-mail contact: info@openlb.net
5 * The most recent release of OpenLB can be downloaded at
6 * <http://www.openlb.net/>
7 *
8 * This program is free software; you can redistribute it and/or
9 * modify it under the terms of the GNU General Public License
10 * as published by the Free Software Foundation; either version 2
11 * of the License, or (at your option) any later version.
12 *
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public
19 * License along with this program; if not, write to the Free
20 * Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
21 * Boston, MA 02110-1301, USA.
22 */
23
24#ifndef VTU_SURFACE_WRITER_H
25#define VTU_SURFACE_WRITER_H
26
29#include "io/base64.hh"
30#include "io/stlReader.h"
32
33namespace olb {
34
35template <typename T>
37protected:
38 std::string _name;
39 std::vector<SuperF3D<T, T>*> _functors;
40 std::vector<AnalyticalF3D<T, T>*> _functorsA;
41 std::vector<Vector<T, 3>> pos;
42 std::vector<std::vector<int>> connectivity;
45 bool _binary;
46
48
49 // performes <VTKFile ...> and <Collection>
50 void preamblePVD(const std::string& fullNamePVD);
51 // performes </Collection> and </VTKFile>
52 void closePVD(const std::string& fullNamePVD);
54 void preambleVTU(const std::string& fullName, int num);
56 void closeVTU(const std::string& fullNamePiece);
58 void dataPVD(int iT, int i, const std::string& fullNamePVD, const std::string& namePiece);
60 void dataPVDmaster(int iT, int i, const std::string& fullNamePVDMaster, const std::string& namePiece);
61 void dataArray(const std::string& fullName);
63 void writeFunctor(const std::string& fullName, std::ofstream& fout, SuperF3D<T, T>& f);
65 void writeAnalyticalFunctor(const std::string& fullName, std::ofstream& fout, AnalyticalF3D<T, T>& f);
67 void writePosition(std::ofstream& fout);
68 static bool areVectorsEqual(const std::vector<T>& v1, const std::vector<T>& v2, T epsilon = 1e-6);
70 std::vector<std::vector<T>> removeDuplicates(const std::vector<std::vector<T>>& points);
71 std::vector<std::vector<int>> mapTrianglesToIndices(STLreader<T>& stlReader);
72
73public:
74 VTUsurfaceWriter(std::string name, CuboidDecomposition3D<T>& cuboidDecomposition, LoadBalancer<T>& loadBalancer,
75 bool binary = false);
76
77 void createMasterFile();
78 void write(std::size_t iT);
79
81 void addPoint(Vector<T, 3>& new_position);
83 void addPoints(std::vector<Vector<T, 3>>& new_positions);
85 void addConnectivity(std::vector<std::vector<int>>& connections);
86
88 void addSTL(STLreader<T>& stlReader);
89
94
95};
96
97} // namespace olb
98
99#endif
AnalyticalF are applications from DD to XD, where X is set by the constructor.
Base class for all LoadBalancer.
Definition vtiWriter.h:42
class for marking output with some text
represents all functors that operate on a SuperStructure<T,3> in general
Definition aliases.h:184
void addPoints(std::vector< Vector< T, 3 > > &new_positions)
Add multiple locations to output set.
void addPoint(Vector< T, 3 > &new_position)
Add single location to output set.
void writeAnalyticalFunctor(const std::string &fullName, std::ofstream &fout, AnalyticalF3D< T, T > &f)
writes functors.
void addConnectivity(std::vector< std::vector< int > > &connections)
Define connectivity between previously-added output points, rendering them a surface.
void dataPVD(int iT, int i, const std::string &fullNamePVD, const std::string &namePiece)
performes <DataSet timestep= ... file=namePiece >
static bool areVectorsEqual(const std::vector< T > &v1, const std::vector< T > &v2, T epsilon=1e-6)
void closeVTU(const std::string &fullNamePiece)
performes </ImageData> and </VTKFile>
VTUsurfaceWriter(std::string name, CuboidDecomposition3D< T > &cuboidDecomposition, LoadBalancer< T > &loadBalancer, bool binary=false)
std::vector< Vector< T, 3 > > pos
void addFunctor(AnalyticalF3D< T, T > &f)
Schedule analytical functor for output.
CuboidDecomposition3D< T > & _cuboidDecomposition
std::vector< AnalyticalF3D< T, T > * > _functorsA
std::vector< SuperF3D< T, T > * > _functors
void write(std::size_t iT)
writes vtu files and assigns values
void closePVD(const std::string &fullNamePVD)
void writePosition(std::ofstream &fout)
writes coordinates of points
void dataPVDmaster(int iT, int i, const std::string &fullNamePVDMaster, const std::string &namePiece)
performes <DataSet timestep= ... file=namePiece >
std::vector< std::vector< T > > removeDuplicates(const std::vector< std::vector< T > > &points)
removes duplicates from points, needed for STL files
void preambleVTU(const std::string &fullName, int num)
performes <VTKFile ...>, <ImageData ...> and <PieceExtent ...>
void addSTL(STLreader< T > &stlReader)
Add full STL surface as output surface.
void writeFunctor(const std::string &fullName, std::ofstream &fout, SuperF3D< T, T > &f)
interpolates and writes functors stored at functors
std::vector< std::vector< int > > connectivity
void preamblePVD(const std::string &fullNamePVD)
std::vector< std::vector< int > > mapTrianglesToIndices(STLreader< T > &stlReader)
void dataArray(const std::string &fullName)
LoadBalancer< T > & _loadBalancer
Plain old scalar vector.
Wrapper functions that simplify the use of MPI.
Top level namespace for all of OpenLB.
Input in STL format – header file.