writer.hh 6.3 KB
Newer Older
dedner's avatar
dedner committed
1
2
3
4
5
6
7
8
9
// -*- tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:

#ifndef DUNE_PYTHON_VTK_WRITER_HH
#define DUNE_PYTHON_VTK_WRITER_HH

#include <dune/vtk/vtkwriter.hh>

#include <dune/python/pybind11/pybind11.h>
10
#include <dune/python/pybind11/stl.h>
dedner's avatar
dedner committed
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27

namespace Dune
{

  namespace Vtk
  {

    // registerVTKWriter
    // -----------------

    template< class Writer, class... options >
    inline static void registerVtkWriter ( pybind11::handle scope,
                                           pybind11::class_< Writer, options... > cls )
    {
      using GridView = typename Writer::GridView;
      using VirtualizedGF = Dune::Vtk::Function<GridView>;

28
29
30
31
32
33
34
35
36
37
38
      cls.def( pybind11::init( [] ( GridView &grid,
               Vtk::FormatTypes format,
               Vtk::DataTypes datatype,
               Vtk::DataTypes headertype,
               pybind11::kwargs kwargs) {
        return new Writer( grid, format, datatype, headertype ); 
        }), pybind11::arg("grid"),
            pybind11::arg("format") = Vtk::FormatTypes::BINARY,
            pybind11::arg("datatype") = Vtk::DataTypes::FLOAT32,
            pybind11::arg("headertype") = Vtk::DataTypes::UINT32,
            pybind11::keep_alive< 1, 2 >());
dedner's avatar
dedner committed
39
40

      cls.def( "write",
41
          [] ( Writer &writer, const std::string &name ) {
dedner's avatar
dedner committed
42
43
44
45
46
47
48
49
50
51
            writer.write( name );
          },
          pybind11::arg("name") );
      cls.def( "write",
          [] ( Writer &writer, const std::string &name, int number ) {
            std::stringstream s; s << name << std::setw(5) << std::setfill('0') << number;
            writer.write( s.str());
          },
          pybind11::arg("name"),
          pybind11::arg("number") );
52
53
54
55
56

      cls.def( "addPointData",
          [] ( Writer &writer, VirtualizedGF &f,
               RangeTypes range, DataTypes data
             ) {
57
58
            f.setRangeType(range);
            f.setDataType(data);
59
60
            writer.addPointData(f);
          },
61
          pybind11::keep_alive< 1, 2 >(),
62
63
64
65
66
67
68
69
70
71
          pybind11::arg("f"),
          pybind11::arg("range")=RangeTypes::AUTO,
          pybind11::arg("data")=DataTypes::FLOAT32
        );
      cls.def( "addPointData",
          [] ( Writer &writer, VirtualizedGF &f,
               std::string &name,
               RangeTypes range, DataTypes data
             ) {
            f.setName(name);
72
73
            f.setRangeType(range);
            f.setDataType(data);
74
75
            writer.addPointData(f);
          },
76
          pybind11::keep_alive< 1, 2 >(),
77
78
79
80
81
82
83
84
85
86
87
          pybind11::arg("f"), pybind11::arg("name"),
          pybind11::arg("range")=RangeTypes::AUTO,
          pybind11::arg("data")=DataTypes::FLOAT32
        );
      cls.def( "addPointData",
          [] ( Writer &writer, VirtualizedGF &f,
               std::string &name, 
               std::vector<int> &components,
               RangeTypes range, DataTypes data
             ) {
            f.setName(name);
88
89
            f.setRangeType(range);
            f.setDataType(data);
90
91
92
            f.setComponents(components);
            writer.addPointData(f);
          },
93
          pybind11::keep_alive< 1, 2 >(),
94
95
96
97
98
99
100
101
102
103
104
          pybind11::arg("f"), pybind11::arg("name"),
          pybind11::arg("components"),
          pybind11::arg("range")=RangeTypes::AUTO,
          pybind11::arg("data")=DataTypes::FLOAT32
        );
      cls.def( "addPointData",
          [] ( Writer &writer, VirtualizedGF &f,
               FieldInfo &info) {
            f.setFieldInfo(info);
            writer.addPointData(f);
          },
105
          pybind11::keep_alive< 1, 2 >(),
106
          pybind11::arg("f"), pybind11::arg("info") );
dedner's avatar
dedner committed
107
      cls.def( "addPointData",
108
109
110
111
          [] ( Writer &writer, VirtualizedGF &f,
               std::vector<int> &components, FieldInfo &info ) {
            f.setFieldInfo(info);
            f.setComponents(components);
dedner's avatar
dedner committed
112
113
            writer.addPointData(f);
          },
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
          pybind11::keep_alive< 1, 2 >(),
          pybind11::arg("f"), pybind11::arg("components"), pybind11::arg("info") );

      cls.def( "addCellData",
          [] ( Writer &writer, VirtualizedGF &f,
               RangeTypes range, DataTypes data
             ) {
            f.setRangeType(range);
            f.setDataType(data);
            writer.addCellData(f);
          },
          pybind11::keep_alive< 1, 2 >(),
          pybind11::arg("f"),
          pybind11::arg("range")=RangeTypes::AUTO,
          pybind11::arg("data")=DataTypes::FLOAT32
        );
      cls.def( "addCellData",
          [] ( Writer &writer, VirtualizedGF &f,
               std::string &name,
               RangeTypes range, DataTypes data
             ) {
            f.setName(name);
            f.setRangeType(range);
            f.setDataType(data);
            writer.addCellData(f);
          },
          pybind11::keep_alive< 1, 2 >(),
          pybind11::arg("f"), pybind11::arg("name"),
          pybind11::arg("range")=RangeTypes::AUTO,
          pybind11::arg("data")=DataTypes::FLOAT32
        );
      cls.def( "addCellData",
          [] ( Writer &writer, VirtualizedGF &f,
               std::string &name, 
               std::vector<int> &components,
               RangeTypes range, DataTypes data
             ) {
            f.setName(name);
            f.setRangeType(range);
            f.setDataType(data);
            f.setComponents(components);
            writer.addCellData(f);
          },
          pybind11::keep_alive< 1, 2 >(),
          pybind11::arg("f"), pybind11::arg("name"),
          pybind11::arg("components"),
          pybind11::arg("range")=RangeTypes::AUTO,
          pybind11::arg("data")=DataTypes::FLOAT32
        );
      cls.def( "addCellData",
          [] ( Writer &writer, VirtualizedGF &f,
               FieldInfo &info) {
            f.setFieldInfo(info);
            writer.addCellData(f);
          },
          pybind11::keep_alive< 1, 2 >(),
          pybind11::arg("f"), pybind11::arg("info") );
      cls.def( "addCellData",
          [] ( Writer &writer, VirtualizedGF &f,
               std::vector<int> &components, FieldInfo &info ) {
            f.setFieldInfo(info);
            f.setComponents(components);
            writer.addCellData(f);
          },
          pybind11::keep_alive< 1, 2 >(),
179
          pybind11::arg("f"), pybind11::arg("components"), pybind11::arg("info") );
dedner's avatar
dedner committed
180
181
    }

182
  } // namespace Vtk
dedner's avatar
dedner committed
183
184
185

} // namespace Dune

186
#endif // DUNE_PYTHON_VTK_WRITER_HH