Merge lp:~nschloe/dolfin/exodus into lp:~fenics-core/dolfin/trunk
- exodus
- Merge into trunk
Status: | Merged |
---|---|
Merged at revision: | 7521 |
Proposed branch: | lp:~nschloe/dolfin/exodus |
Merge into: | lp:~fenics-core/dolfin/trunk |
Diff against target: |
598 lines (+536/-0) 7 files modified
dolfin/CMakeLists.txt (+15/-0) dolfin/common/defines.cpp (+9/-0) dolfin/common/defines.h (+3/-0) dolfin/io/ExodusFile.cpp (+335/-0) dolfin/io/ExodusFile.h (+75/-0) dolfin/io/File.cpp (+7/-0) test/unit/io/python/Exodus.py (+92/-0) |
To merge this branch: | bzr merge lp:~nschloe/dolfin/exodus |
Related bugs: |
Reviewer | Review Type | Date Requested | Status |
---|---|---|---|
Garth Wells | Pending | ||
Review via email: mp+151023@code.launchpad.net |
Commit message
Description of the change
This is essentially an update of a former attempt at getting proper Exodus support into Dolfin, with the exception that it behaves gracefully if the present VTK installation doesn't have Exodus support.
- 7482. By Nico Schlömer
-
actually add defintion -DHAS_VTK_EXODUS
- 7483. By Nico Schlömer
-
nest HAS_VTK and HAS_VTK_EXODUS
- 7484. By Nico Schlömer
-
Unit* -> Unit*Mesh in Exodus tests
Anders Logg (logg) wrote : | # |
- 7485. By Nico Schlömer
-
merge with HEAD
Nico Schlömer (nschloe) wrote : | # |
It's merged now, still looking good.
Anders Logg (logg) wrote : | # |
On Thu, Mar 07, 2013 at 10:58:22PM -0000, Nico Schlömer wrote:
> It's merged now, still looking good.
ok, I'll try looking at it tomorrow. It might now make 1.2 but there
will be new releases.
--
Anders
Nico Schlömer (nschloe) wrote : | # |
Shorter release cycles ftw. https:/
Anders Logg (logg) wrote : | # |
Doesn't seem to work for me. Did you forget to merge in File.cpp?
Error: Unable to open file.
*** Reason: Unknown file type (".e") for file "mesh.e".
*** Where: This error was encountered inside File.cpp.
*** Process: 0
- 7486. By Nico Schlömer
-
merge with HEAD
- 7487. By Nico Schlömer
-
added CMake status message if vtkExodusWriter isn't found
Nico Schlömer (nschloe) wrote : | # |
Exodus file writing is only activated if the system VTK has Exodus support. There's some CMake logic to check for it, and I just added a status message if it's not found (which may be the case on your system). With the default Ubuntu VTK, for example, it does work.
If you want to try it out again, you should either see a configure time status message about missing Exodus or it should work.
Anders Logg (logg) wrote : | # |
On Mon, Mar 11, 2013 at 12:52:24PM -0000, Nico Schlömer wrote:
> Exodus file writing is only activated if the system VTK has Exodus support. There's some CMake logic to check for it, and I just added a status message if it's not found (which may be the case on your system). With the default Ubuntu VTK, for example, it does work.
> If you want to try it out again, you should either see a configure time status message about missing Exodus or it should work.
The unit test should work regardless of whether VTK has Exodus
support.
I suggest adding:
1. a runtime error (dolfin_error) if VTK Exodus is missing (in
File.cpp) instead of just missing file format
2. a runtime check in the unit test so that it is skipped if the VTK
Exodus format is missing
--
Anders
- 7488. By Nico Schlömer
-
added dolfin:
:has_exodus( ) and include Exodus unit tests conditionally on this
Nico Schlömer (nschloe) wrote : | # |
Ah, the unit test, of course.
I now went the same way that HDF5 and XDMF support is handled: defines.cpp now features a has_exodus() helper function, and the unit tests are added only if it returns `true`.
Anders Logg (logg) wrote : | # |
On Mon, Mar 11, 2013 at 01:16:17PM -0000, Nico Schlömer wrote:
> Ah, the unit test, of course.
> I now went the same way that HDF5 and XDMF support is handled: defines.cpp now features a has_exodus() helper function, and the unit tests are added only if it returns `true`.
Thanks, looks good now.
--
Anders
Preview Diff
1 | === modified file 'dolfin/CMakeLists.txt' |
2 | --- dolfin/CMakeLists.txt 2013-03-02 20:24:12 +0000 |
3 | +++ dolfin/CMakeLists.txt 2013-03-11 13:15:27 +0000 |
4 | @@ -250,6 +250,21 @@ |
5 | endif() |
6 | endif() |
7 | |
8 | + # Check if the VTK-Exodus library is found. |
9 | + # Normally, it should be included in libvtkParallel.so (check out the output |
10 | + # of |
11 | + # $ nm -D /usr/lib/libvtkParallel.so | grep ExodusIIWriter |
12 | + # but on some systems, this doesn't seem to be the case. Check if the |
13 | + # appropriate header exists and store this information in HAS_VTK_EXODUS. |
14 | + include(CheckIncludeFileCXX) |
15 | + set(CMAKE_REQUIRED_INCLUDES "${VTK_INCLUDE_DIRS}") |
16 | + CHECK_INCLUDE_FILE_CXX(vtkExodusIIWriter.h HAS_VTK_EXODUS) |
17 | + if(HAS_VTK_EXODUS) |
18 | + list(APPEND DOLFIN_CXX_DEFINITIONS "-DHAS_VTK_EXODUS") |
19 | + else() |
20 | + message(STATUS "vtkExodusIIWriter not found, disabling Exodus file writing") |
21 | + endif() |
22 | + |
23 | list(APPEND DOLFIN_TARGET_LINK_LIBRARIES ${DOLFIN_VTK_LIBRARIES}) |
24 | endif() |
25 | |
26 | |
27 | === modified file 'dolfin/common/defines.cpp' |
28 | --- dolfin/common/defines.cpp 2012-09-14 08:39:20 +0000 |
29 | +++ dolfin/common/defines.cpp 2013-03-11 13:15:27 +0000 |
30 | @@ -122,3 +122,12 @@ |
31 | #endif |
32 | } |
33 | //------------------------------------------------------------------------- |
34 | +bool dolfin::has_exodus() |
35 | +{ |
36 | +#ifdef HAS_VTK_EXODUS |
37 | + return true; |
38 | +#else |
39 | + return false; |
40 | +#endif |
41 | +} |
42 | +//------------------------------------------------------------------------- |
43 | |
44 | === modified file 'dolfin/common/defines.h' |
45 | --- dolfin/common/defines.h 2012-09-14 08:39:20 +0000 |
46 | +++ dolfin/common/defines.h 2013-03-11 13:15:27 +0000 |
47 | @@ -61,6 +61,9 @@ |
48 | // Return true if DOLFIN is compiled with HDF5 |
49 | bool has_hdf5(); |
50 | |
51 | + // Return true if DOLFIN is compiled with Exodus |
52 | + bool has_exodus(); |
53 | + |
54 | } |
55 | |
56 | #endif |
57 | |
58 | === added file 'dolfin/io/ExodusFile.cpp' |
59 | --- dolfin/io/ExodusFile.cpp 1970-01-01 00:00:00 +0000 |
60 | +++ dolfin/io/ExodusFile.cpp 2013-03-11 13:15:27 +0000 |
61 | @@ -0,0 +1,335 @@ |
62 | +// Copyright (C) 2013 Nico Schlömer |
63 | +// |
64 | +// This file is part of DOLFIN. |
65 | +// |
66 | +// DOLFIN is free software: you can redistribute it and/or modify |
67 | +// it under the terms of the GNU Lesser General Public License as published by |
68 | +// the Free Software Foundation, either version 3 of the License, or |
69 | +// (at your option) any later version. |
70 | +// |
71 | +// DOLFIN is distributed in the hope that it will be useful, |
72 | +// but WITHOUT ANY WARRANTY; without even the implied warranty of |
73 | +// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
74 | +// GNU Lesser General Public License for more details. |
75 | +// |
76 | +// You should have received a copy of the GNU Lesser General Public License |
77 | +// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>. |
78 | +// |
79 | +// First added: 2013-02-27 |
80 | + |
81 | +#ifdef HAS_VTK |
82 | +#ifdef HAS_VTK_EXODUS |
83 | + |
84 | +#include <dolfin/fem/GenericDofMap.h> |
85 | +#include <dolfin/function/Function.h> |
86 | +#include <dolfin/function/FunctionSpace.h> |
87 | +#include <dolfin/la/GenericVector.h> |
88 | +#include <dolfin/mesh/Mesh.h> |
89 | +#include <dolfin/mesh/MeshFunction.h> |
90 | +#include "ExodusFile.h" |
91 | + |
92 | +#include <vtkUnsignedIntArray.h> |
93 | +#include <vtkIntArray.h> |
94 | +#include <vtkDoubleArray.h> |
95 | +#include <vtkCellType.h> |
96 | +#include <vtkPointData.h> |
97 | +#include <vtkCellData.h> |
98 | +#include <vtkCellArray.h> |
99 | +#include <vtkIdTypeArray.h> |
100 | +#include <vtkUnstructuredGrid.h> |
101 | +#include <vtkExodusIIWriter.h> |
102 | + |
103 | +using namespace dolfin; |
104 | + |
105 | +//---------------------------------------------------------------------------- |
106 | +ExodusFile::ExodusFile(const std::string filename) |
107 | + : GenericFile(filename, "Exodus") |
108 | +{ |
109 | +} |
110 | +//---------------------------------------------------------------------------- |
111 | +ExodusFile::~ExodusFile() |
112 | +{ |
113 | + // Do nothing |
114 | +} |
115 | +//---------------------------------------------------------------------------- |
116 | +void ExodusFile::operator<<(const Mesh& mesh) |
117 | +{ |
118 | + perform_write(create_vtk_mesh(mesh)); |
119 | + log(TRACE, "Saved mesh %s (%s) to file %s in Exodus format.", |
120 | + mesh.name().c_str(), mesh.label().c_str(), _filename.c_str()); |
121 | + return; |
122 | +} |
123 | +//---------------------------------------------------------------------------- |
124 | +void ExodusFile::operator<<(const MeshFunction<unsigned int>& meshfunction) |
125 | +{ |
126 | + const Mesh& mesh = meshfunction.mesh(); |
127 | + const uint cell_dim = meshfunction.dim(); |
128 | + |
129 | + // Throw error for MeshFunctions on vertices for interval elements |
130 | + if (mesh.topology().dim() == 1 && cell_dim == 0) |
131 | + { |
132 | + dolfin_error("ExodusFile.cpp", |
133 | + "write mesh function to Exodus file", |
134 | + "Exodus output of mesh functions on interval facets is not supported"); |
135 | + } |
136 | + |
137 | + if (cell_dim != mesh.topology().dim() && cell_dim != mesh.topology().dim() - 1) |
138 | + { |
139 | + dolfin_error("ExodusFile.cpp", |
140 | + "write mesh function to Exodus file", |
141 | + "Exodus output of mesh functions is implemented for cell- and facet-based functions only"); |
142 | + } |
143 | + |
144 | + // Create Exodus mesh. |
145 | + vtkSmartPointer<vtkUnstructuredGrid> vtk_mesh = create_vtk_mesh(mesh); |
146 | + |
147 | + // Add cell data. |
148 | + const int dim = meshfunction.dim(); |
149 | + const int numCells = mesh.num_cells(); |
150 | + vtkSmartPointer<vtkUnsignedIntArray> cellData = |
151 | + vtkSmartPointer<vtkUnsignedIntArray>::New(); |
152 | + cellData->SetNumberOfComponents(dim); |
153 | + cellData->SetArray(const_cast<unsigned int*>(meshfunction.values()), dim*numCells, 1); |
154 | + cellData->SetName(meshfunction.name().c_str()); |
155 | + vtk_mesh->GetCellData()->AddArray(cellData); |
156 | + |
157 | + // Write out. |
158 | + perform_write(vtk_mesh); |
159 | + |
160 | + log(TRACE, "Saved mesh function %s (%s) to file %s in Exodus format.", |
161 | + mesh.name().c_str(), mesh.label().c_str(), _filename.c_str()); |
162 | +} |
163 | +//---------------------------------------------------------------------------- |
164 | +void ExodusFile::operator<<(const MeshFunction<int>& meshfunction) |
165 | +{ |
166 | + const Mesh& mesh = meshfunction.mesh(); |
167 | + const uint cell_dim = meshfunction.dim(); |
168 | + |
169 | + // Throw error for MeshFunctions on vertices for interval elements |
170 | + if (mesh.topology().dim() == 1 && cell_dim == 0) |
171 | + { |
172 | + dolfin_error("ExodusFile.cpp", |
173 | + "write mesh function to Exodus file", |
174 | + "Exodus output of mesh functions on interval facets is not supported"); |
175 | + } |
176 | + |
177 | + if (cell_dim != mesh.topology().dim() && cell_dim != mesh.topology().dim() - 1) |
178 | + { |
179 | + dolfin_error("ExodusFile.cpp", |
180 | + "write mesh function to Exodus file", |
181 | + "Exodus output of mesh functions is implemented for cell- and facet-based functions only"); |
182 | + } |
183 | + |
184 | + // Create Exodus mesh. |
185 | + vtkSmartPointer<vtkUnstructuredGrid> vtk_mesh = create_vtk_mesh(mesh); |
186 | + |
187 | + // Add cell data. |
188 | + const int dim = meshfunction.dim(); |
189 | + const int numCells = mesh.num_cells(); |
190 | + vtkSmartPointer<vtkIntArray> cellData = |
191 | + vtkSmartPointer<vtkIntArray>::New(); |
192 | + cellData->SetNumberOfComponents(dim); |
193 | + cellData->SetArray(const_cast<int*>(meshfunction.values()), dim*numCells, 1); |
194 | + cellData->SetName(meshfunction.name().c_str()); |
195 | + vtk_mesh->GetCellData()->AddArray(cellData); |
196 | + |
197 | + // Write out. |
198 | + perform_write(vtk_mesh); |
199 | + |
200 | + log(TRACE, "Saved mesh function %s (%s) to file %s in Exodus format.", |
201 | + mesh.name().c_str(), mesh.label().c_str(), _filename.c_str()); |
202 | +} |
203 | +//---------------------------------------------------------------------------- |
204 | +void ExodusFile::operator<<(const MeshFunction<double>& meshfunction) |
205 | +{ |
206 | + const Mesh& mesh = meshfunction.mesh(); |
207 | + const uint cell_dim = meshfunction.dim(); |
208 | + |
209 | + // Throw error for MeshFunctions on vertices for interval elements |
210 | + if (mesh.topology().dim() == 1 && cell_dim == 0) |
211 | + { |
212 | + dolfin_error("ExodusFile.cpp", |
213 | + "write mesh function to Exodus file", |
214 | + "Exodus output of mesh functions on interval facets is not supported"); |
215 | + } |
216 | + |
217 | + if (cell_dim != mesh.topology().dim() && cell_dim != mesh.topology().dim() - 1) |
218 | + { |
219 | + dolfin_error("ExodusFile.cpp", |
220 | + "write mesh function to Exodus file", |
221 | + "Exodus output of mesh functions is implemented for cell- and facet-based functions only"); |
222 | + } |
223 | + |
224 | + // Create Exodus mesh. |
225 | + vtkSmartPointer<vtkUnstructuredGrid> vtk_mesh = create_vtk_mesh(mesh); |
226 | + |
227 | + // Add cell data. |
228 | + const int dim = meshfunction.dim(); |
229 | + const int numCells = mesh.num_cells(); |
230 | + vtkSmartPointer<vtkDoubleArray> cellData = |
231 | + vtkSmartPointer<vtkDoubleArray>::New(); |
232 | + cellData->SetNumberOfComponents(dim); |
233 | + cellData->SetArray(const_cast<double*>(meshfunction.values()), dim*numCells, 1); |
234 | + cellData->SetName(meshfunction.name().c_str()); |
235 | + vtk_mesh->GetCellData()->AddArray(cellData); |
236 | + |
237 | + // Write out. |
238 | + perform_write(vtk_mesh); |
239 | + |
240 | + log(TRACE, "Saved mesh function %s (%s) to file %s in Exodus format.", |
241 | + mesh.name().c_str(), mesh.label().c_str(), _filename.c_str()); |
242 | +} |
243 | +//---------------------------------------------------------------------------- |
244 | +void ExodusFile::operator<<(const Function& u) |
245 | +{ |
246 | + u.update(); |
247 | + write_function(u, counter); |
248 | +} |
249 | +//---------------------------------------------------------------------------- |
250 | +void ExodusFile::operator<<(const std::pair<const Function*, double> u) |
251 | +{ |
252 | + dolfin_assert(u.first); |
253 | + u.first->update(); |
254 | + write_function(*(u.first), u.second); |
255 | +} |
256 | +//---------------------------------------------------------------------------- |
257 | +void ExodusFile::write_function(const Function& u, double time) const |
258 | +{ |
259 | + // Write results |
260 | + // Get rank of Function |
261 | + const uint rank = u.value_rank(); |
262 | + if (rank > 2) |
263 | + { |
264 | + dolfin_error("ExodusFile.cpp", |
265 | + "write data to Exodus file", |
266 | + "Only scalar, vector and tensor functions can be saved in Exodus format"); |
267 | + } |
268 | + |
269 | + // Get number of components |
270 | + const uint dim = u.value_size(); |
271 | + |
272 | + // Check that function type can be handled, cf. |
273 | + // http://www.vtk.org/Bug/view.php?id=13508. |
274 | + if (dim > 6) |
275 | + { |
276 | + dolfin_error("ExodusFile.cpp", |
277 | + "write data to Exodus file", |
278 | + "Can't handle more than 6 components"); |
279 | + } |
280 | + |
281 | + // Test for cell-based element type |
282 | + dolfin_assert(u.function_space()->mesh()); |
283 | + const Mesh& mesh = *u.function_space()->mesh(); |
284 | + |
285 | + vtkSmartPointer<vtkUnstructuredGrid> vtk_mesh = |
286 | + create_vtk_mesh(mesh); |
287 | + |
288 | + uint cell_based_dim = 1; |
289 | + for (uint i = 0; i < rank; i++) |
290 | + cell_based_dim *= mesh.topology().dim(); |
291 | + |
292 | + dolfin_assert(u.function_space()->dofmap()); |
293 | + const GenericDofMap& dofmap = *u.function_space()->dofmap(); |
294 | + // Define the vector that holds the values outside the if |
295 | + // to make sure it doesn't get destroyed before the data |
296 | + // is written out to a file |
297 | + std::vector<double> values; |
298 | + if (dofmap.max_cell_dimension() == cell_based_dim) |
299 | + { |
300 | + // Extract DOFs from u. |
301 | + const uint num_cells = mesh.num_cells(); |
302 | + const uint size = num_cells*dim; |
303 | + std::vector<int> dof_set; |
304 | + for (CellIterator cell(mesh); !cell.end(); ++cell) |
305 | + { |
306 | + const std::vector<int>& dofs = dofmap.cell_dofs(cell->index()); |
307 | + for(uint i = 0; i < dofmap.cell_dimension(cell->index()); ++i) |
308 | + dof_set.push_back(dofs[i]); |
309 | + } |
310 | + // Get values |
311 | + values.resize(dof_set.size()); |
312 | + dolfin_assert(u.vector()); |
313 | + u.vector()->get_local(&values[0], dof_set.size(), &dof_set[0]); |
314 | + |
315 | + // Set the cell array. |
316 | + vtkSmartPointer<vtkDoubleArray> cellData = |
317 | + vtkSmartPointer<vtkDoubleArray>::New(); |
318 | + cellData->SetNumberOfComponents(dim); |
319 | + cellData->SetArray(&values[0], dof_set.size(), 1); |
320 | + cellData->SetName(u.name().c_str()); |
321 | + vtk_mesh->GetCellData()->AddArray(cellData); |
322 | + } |
323 | + else |
324 | + { |
325 | + // Extract point values. |
326 | + const uint num_vertices = mesh.num_vertices(); |
327 | + const uint size = num_vertices*dim; |
328 | + values.resize(size); |
329 | + u.compute_vertex_values(values, mesh); |
330 | + // Set the point array. |
331 | + vtkSmartPointer<vtkDoubleArray> pointData = |
332 | + vtkSmartPointer<vtkDoubleArray>::New(); |
333 | + pointData->SetNumberOfComponents(dim); |
334 | + pointData->SetArray(&values[0], size, 1); |
335 | + pointData->SetName(u.name().c_str()); |
336 | + vtk_mesh->GetPointData()->AddArray(pointData); |
337 | + } |
338 | + |
339 | + // Actually write out the data. |
340 | + perform_write(vtk_mesh); |
341 | + |
342 | + log(TRACE, "Saved function %s (%s) to file %s in Exodus format.", |
343 | + u.name().c_str(), u.label().c_str(), _filename.c_str()); |
344 | +} |
345 | +//---------------------------------------------------------------------------- |
346 | +vtkSmartPointer<vtkUnstructuredGrid> ExodusFile::create_vtk_mesh(const Mesh& mesh) const |
347 | +{ |
348 | + // Build Exodus unstructured grid object. |
349 | + vtkSmartPointer<vtkUnstructuredGrid> unstructuredGrid = |
350 | + vtkSmartPointer<vtkUnstructuredGrid>::New(); |
351 | + // Set the points. |
352 | + const int numPoints = mesh.num_vertices(); |
353 | + vtkSmartPointer<vtkDoubleArray> pointData = |
354 | + vtkSmartPointer<vtkDoubleArray>::New(); |
355 | + pointData->SetNumberOfComponents(3); |
356 | + pointData->SetArray(const_cast<double*>(&mesh.coordinates()[0]), 3*numPoints, 1); |
357 | + vtkSmartPointer<vtkPoints> points = |
358 | + vtkSmartPointer<vtkPoints>::New(); |
359 | + points->SetData(pointData); |
360 | + unstructuredGrid->SetPoints(points); |
361 | + // Set cells. Those need to be copied over since the |
362 | + // default Dolfin node ID data type is dolfin::uint |
363 | + // (typically unsigned int), and the node ID of Exodus |
364 | + // is vtkIdType (typically long long int). |
365 | + const int numCells = mesh.num_cells(); |
366 | + const std::vector<unsigned int> cells = mesh.cells(); |
367 | + vtkSmartPointer<vtkCellArray> cellData = |
368 | + vtkSmartPointer<vtkCellArray>::New(); |
369 | + vtkIdType tmp[4]; |
370 | + for (int k=0; k<numCells; k++) |
371 | + { |
372 | + for (int i=0; i<4; i++) |
373 | + tmp[i] = cells[4*k+i]; |
374 | + cellData->InsertNextCell(4, tmp); |
375 | + } |
376 | + unstructuredGrid->SetCells(VTK_TETRA, cellData); |
377 | + |
378 | + return unstructuredGrid; |
379 | +} |
380 | +//---------------------------------------------------------------------------- |
381 | +void ExodusFile::perform_write(const vtkSmartPointer<vtkUnstructuredGrid> & vtk_mesh) const |
382 | +{ |
383 | + // Instantiate writer. |
384 | + vtkSmartPointer<vtkExodusIIWriter> writer = |
385 | + vtkSmartPointer<vtkExodusIIWriter>::New(); |
386 | + |
387 | + // Write out to file. |
388 | + writer->SetFileName(_filename.c_str()); |
389 | + writer->SetInput(vtk_mesh); |
390 | + writer->Write(); |
391 | + |
392 | + return; |
393 | +} |
394 | +//---------------------------------------------------------------------------- |
395 | +#endif // HAS_VTK_EXODUS |
396 | +#endif // HAS_VTK |
397 | |
398 | === added file 'dolfin/io/ExodusFile.h' |
399 | --- dolfin/io/ExodusFile.h 1970-01-01 00:00:00 +0000 |
400 | +++ dolfin/io/ExodusFile.h 2013-03-11 13:15:27 +0000 |
401 | @@ -0,0 +1,75 @@ |
402 | +// Copyright (C) 2013 Nico Schlömer |
403 | +// |
404 | +// This file is part of DOLFIN. |
405 | +// |
406 | +// DOLFIN is free software: you can redistribute it and/or modify |
407 | +// it under the terms of the GNU Lesser General Public License as published by |
408 | +// the Free Software Foundation, either version 3 of the License, or |
409 | +// (at your option) any later version. |
410 | +// |
411 | +// DOLFIN is distributed in the hope that it will be useful, |
412 | +// but WITHOUT ANY WARRANTY; without even the implied warranty of |
413 | +// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
414 | +// GNU Lesser General Public License for more details. |
415 | +// |
416 | +// You should have received a copy of the GNU Lesser General Public License |
417 | +// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>. |
418 | +// |
419 | +// First added: 2013-02-27 |
420 | + |
421 | +#ifndef __EXODUS_FILE_H |
422 | +#define __EXODUS_FILE_H |
423 | + |
424 | +#ifdef HAS_VTK |
425 | +#ifdef HAS_VTK_EXODUS |
426 | + |
427 | +#include "GenericFile.h" |
428 | + |
429 | +#include <vtkSmartPointer.h> |
430 | + |
431 | +// forward declarations |
432 | +class vtkUnstructuredGrid; |
433 | + |
434 | +namespace dolfin |
435 | +{ |
436 | + |
437 | + /// This class supports the output of meshes and functions |
438 | + /// Exodus format for visualistion purposes. It is not suitable to |
439 | + /// checkpointing as it may decimate some data. |
440 | + |
441 | + class ExodusFile : public GenericFile |
442 | + { |
443 | + public: |
444 | + |
445 | + ExodusFile(const std::string filename); |
446 | + ~ExodusFile(); |
447 | + |
448 | + void operator<< (const Mesh& mesh); |
449 | + void operator<< (const MeshFunction<unsigned int>& meshfunction); |
450 | + void operator<< (const MeshFunction<int>& meshfunction); |
451 | + void operator<< (const MeshFunction<double>& meshfunction); |
452 | + void operator<< (const Function& u); |
453 | + void operator<< (const std::pair<const Function*, double> u); |
454 | + |
455 | +private: |
456 | + |
457 | + void write_function(const Function& u, double time) const; |
458 | + |
459 | + vtkSmartPointer<vtkUnstructuredGrid> create_vtk_mesh(const Mesh& mesh) const; |
460 | + |
461 | + void perform_write(const vtkSmartPointer<vtkUnstructuredGrid> & vtk_mesh) const; |
462 | + |
463 | + // File encoding |
464 | + const std::string encoding; |
465 | + std::string encode_string; |
466 | + |
467 | + bool binary; |
468 | + bool compress; |
469 | + |
470 | + }; |
471 | + |
472 | +} |
473 | + |
474 | +#endif // HAS_VTK_EXODUS |
475 | +#endif // HAS_VTK |
476 | +#endif // __EXODUS_FILE_H |
477 | |
478 | === modified file 'dolfin/io/File.cpp' |
479 | --- dolfin/io/File.cpp 2013-03-04 11:44:12 +0000 |
480 | +++ dolfin/io/File.cpp 2013-03-11 13:15:27 +0000 |
481 | @@ -33,6 +33,7 @@ |
482 | #include "BinaryFile.h" |
483 | #include "RAWFile.h" |
484 | #include "VTKFile.h" |
485 | +#include "ExodusFile.h" |
486 | #include "XMLFile.h" |
487 | #include "XYZFile.h" |
488 | #include "XDMFFile.h" |
489 | @@ -73,6 +74,12 @@ |
490 | file.reset(new XMLFile(filename)); |
491 | else if (extension == ".pvd") |
492 | file.reset(new VTKFile(filename, encoding)); |
493 | +#ifdef HAS_VTK |
494 | +#ifdef HAS_VTK_EXODUS |
495 | + else if (extension == ".e") |
496 | + file.reset(new ExodusFile(filename)); |
497 | +#endif |
498 | +#endif |
499 | else if (extension == ".raw") |
500 | file.reset(new RAWFile(filename)); |
501 | else if (extension == ".xyz") |
502 | |
503 | === added file 'test/unit/io/python/Exodus.py' |
504 | --- test/unit/io/python/Exodus.py 1970-01-01 00:00:00 +0000 |
505 | +++ test/unit/io/python/Exodus.py 2013-03-11 13:15:27 +0000 |
506 | @@ -0,0 +1,92 @@ |
507 | +# -*- coding: utf-8 -*- |
508 | +"""Unit tests for the Exodus io library""" |
509 | + |
510 | +# Copyright (C) 2013 Nico Schlömer |
511 | +# |
512 | +# This file is part of DOLFIN. |
513 | +# |
514 | +# DOLFIN is free software: you can redistribute it and/or modify |
515 | +# it under the terms of the GNU Lesser General Public License as published by |
516 | +# the Free Software Foundation, either version 3 of the License, or |
517 | +# (at your option) any later version. |
518 | +# |
519 | +# DOLFIN is distributed in the hope that it will be useful, |
520 | +# but WITHOUT ANY WARRANTY; without even the implied warranty of |
521 | +# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
522 | +# GNU Lesser General Public License for more details. |
523 | +# |
524 | +# You should have received a copy of the GNU Lesser General Public License |
525 | +# along with DOLFIN. If not, see <http://www.gnu.org/licenses/>. |
526 | +# |
527 | +# First added: 2013-02-27 |
528 | +# Last changed: 2013-03-11 |
529 | + |
530 | +import unittest |
531 | +from dolfin import * |
532 | + |
533 | +if has_exodus(): |
534 | + class Exodus_Mesh_Output(unittest.TestCase): |
535 | + """Test output of Meshes to Exodus files""" |
536 | + |
537 | + def test_save_1d_mesh(self): |
538 | + if MPI.num_processes() == 1: |
539 | + mesh = UnitIntervalMesh(32) |
540 | + File("mesh.e") << mesh |
541 | + |
542 | + def test_save_2d_mesh(self): |
543 | + mesh = UnitSquareMesh(32, 32) |
544 | + File("mesh.e") << mesh |
545 | + |
546 | + def test_save_3d_mesh(self): |
547 | + mesh = UnitCubeMesh(8, 8, 8) |
548 | + File("mesh.e") << mesh |
549 | + |
550 | + |
551 | + class Exodus_Point_Function_Output(unittest.TestCase): |
552 | + """Test output of point-based Functions to Exodus files""" |
553 | + |
554 | + def test_save_1d_scalar(self): |
555 | + if MPI.num_processes() == 1: |
556 | + mesh = UnitIntervalMesh(32) |
557 | + u = Function(FunctionSpace(mesh, "Lagrange", 2)) |
558 | + u.vector()[:] = 1.0 |
559 | + File("u.e") << u |
560 | + |
561 | + def test_save_2d_scalar(self): |
562 | + mesh = UnitSquareMesh(16, 16) |
563 | + u = Function(FunctionSpace(mesh, "Lagrange", 2)) |
564 | + u.vector()[:] = 1.0 |
565 | + File("u.e") << u |
566 | + |
567 | + def test_save_3d_scalar(self): |
568 | + mesh = UnitCubeMesh(8, 8, 8) |
569 | + u = Function(FunctionSpace(mesh, "Lagrange", 2)) |
570 | + u.vector()[:] = 1.0 |
571 | + File("u.e") << u |
572 | + |
573 | + def test_save_2d_vector(self): |
574 | + mesh = UnitSquareMesh(16, 16) |
575 | + u = Function(VectorFunctionSpace(mesh, "Lagrange", 2)) |
576 | + u.vector()[:] = 1.0 |
577 | + File("u.e") << u |
578 | + |
579 | + def test_save_3d_vector(self): |
580 | + mesh = UnitCubeMesh(8, 8, 8) |
581 | + u = Function(VectorFunctionSpace(mesh, "Lagrange", 2)) |
582 | + u.vector()[:] = 1.0 |
583 | + File("u.e") << u |
584 | + |
585 | + def test_save_2d_tensor(self): |
586 | + mesh = UnitSquareMesh(16, 16) |
587 | + u = Function(TensorFunctionSpace(mesh, "Lagrange", 2)) |
588 | + u.vector()[:] = 1.0 |
589 | + File("u.e") << u |
590 | + |
591 | + #def test_save_3d_tensor(self): |
592 | + # mesh = UnitCubeMesh(8, 8, 8) |
593 | + # u = Function(TensorFunctionSpace(mesh, "Lagrange", 2)) |
594 | + # u.vector()[:] = 1.0 |
595 | + # File("u.e") << u |
596 | + |
597 | +if __name__ == "__main__": |
598 | + unittest.main() |
Hi Nico, what's the status of your patch? If you sync it with trunk I can merge it in (if it looks good).