Line data Source code
1 : /*
2 : * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3 : * Author: Andreas Vogel
4 : *
5 : * This file is part of UG4.
6 : *
7 : * UG4 is free software: you can redistribute it and/or modify it under the
8 : * terms of the GNU Lesser General Public License version 3 (as published by the
9 : * Free Software Foundation) with the following additional attribution
10 : * requirements (according to LGPL/GPL v3 §7):
11 : *
12 : * (1) The following notice must be displayed in the Appropriate Legal Notices
13 : * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 : *
15 : * (2) The following notice must be displayed at a prominent place in the
16 : * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 : *
18 : * (3) The following bibliography is recommended for citation and must be
19 : * preserved in all covered files:
20 : * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 : * parallel geometric multigrid solver on hierarchically distributed grids.
22 : * Computing and visualization in science 16, 4 (2013), 151-164"
23 : * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 : * flexible software system for simulating pde based models on high performance
25 : * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 : *
27 : * This program is distributed in the hope that it will be useful,
28 : * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 : * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 : * GNU Lesser General Public License for more details.
31 : */
32 :
33 : #ifndef __H__UG__LIB_DISC__IO__VTKOUTPUT__
34 : #define __H__UG__LIB_DISC__IO__VTKOUTPUT__
35 :
36 : //TODO: deduplicate code in VTUOutput - reduces vulnerability to errors and linker time
37 :
38 : // extern libraries
39 : #include <string>
40 : #include <vector>
41 : #include <map>
42 :
43 : // other ug modules
44 : #include "common/util/string_util.h"
45 : #include "common/util/base64_file_writer.h"
46 : #include "lib_disc/common/function_group.h"
47 : #include "lib_disc/domain.h"
48 : #include "lib_disc/spatial_disc/user_data/user_data.h"
49 :
50 : namespace ug{
51 : // todo this should avoid refactoring all the signatures of VTKOutput, remove it later
52 : typedef Base64FileWriter VTKFileWriter;
53 :
54 : template <typename T>
55 : struct IteratorProvider
56 : {
57 : template <typename TElem>
58 : struct traits
59 : {
60 : typedef typename T::template traits<TElem>::iterator iterator;
61 : typedef typename T::template traits<TElem>::const_iterator const_iterator;
62 : };
63 :
64 : template <typename TElem>
65 0 : static typename traits<TElem>::const_iterator begin(const T& provider, int si)
66 : {
67 0 : if(si < 0) return provider.template begin<TElem>();
68 : else return provider.template begin<TElem>(si);
69 : }
70 :
71 : template <typename TElem>
72 0 : static typename traits<TElem>::const_iterator end(const T& provider, int si)
73 : {
74 0 : if(si < 0) return provider.template end<TElem>();
75 : else return provider.template end<TElem>(si);
76 : }
77 : };
78 :
79 : template <>
80 : struct IteratorProvider<MGSubsetHandler>
81 : {
82 : private:
83 : typedef MGSubsetHandler T;
84 :
85 : public:
86 : template <typename TElem>
87 : struct traits
88 : {
89 : typedef typename T::template traits<TElem>::iterator iterator;
90 : typedef typename T::template traits<TElem>::const_iterator const_iterator;
91 : };
92 :
93 : template <typename TElem>
94 0 : static typename traits<TElem>::const_iterator begin(const T& provider, int si)
95 : {
96 0 : const int lev = provider.num_levels() - 1;
97 0 : if(si < 0) return provider.multi_grid()->template begin<TElem>(lev);
98 0 : else return provider.template begin<TElem>(si, lev);
99 : }
100 :
101 : template <typename TElem>
102 0 : static typename traits<TElem>::const_iterator end(const T& provider, int si)
103 : {
104 0 : const int lev = provider.num_levels() - 1;
105 0 : if(si < 0) return provider.multi_grid()->template end<TElem>(lev);
106 0 : else return provider.template end<TElem>(si, lev);
107 : }
108 : };
109 :
110 :
111 : /// output writer to the VTK file format
112 : /**
113 : * This class can be used to write grid and data associated with the grid to
114 : * the vtk file format.
115 : *
116 : * The produced files are:
117 : *
118 : * 1)) TIME DEPENDENT
119 : * 1a) In case of functions that are only defined for some subsets
120 : * - filename_pXXXX_sXXXX_tXXXX.vtu (piece of grid for proc, subset, timestep)
121 : * - filename_sXXXX_tXXXX.pvtu (piece of grid for subset, timestep)
122 : * - filename_tXXXX.pvd (piece of grid for timestep)
123 : * - filename.pvd (group of timeseries)
124 : *
125 : * 1b) In case that all functions are defined globally
126 : * - filename_p0000_t0000.vtu (piece of grid for proc, timestep)
127 : * - filename_t0000.pvtu (piece of grid for timestep)
128 : *
129 : * 2)) TIME INDEPENDENT
130 : * 2a) In case of functions that are only defined for some subsets
131 : * - filename_pXXXX_sXXXX.vtu (piece of grid for proc, subset)
132 : * - filename_sXXXX.pvtu (piece of grid for subset)
133 : * - filename.pvd (group of subsets)
134 : *
135 : * 2b) In case that all functions are defined globally
136 : * - filename_p0000.vtu (piece of grid for proc)
137 : * - filename.pvtu (group of procs)
138 : *
139 : * In case that the simulation is run in serial (or in parallel environment but
140 : * with only one process), the *.pvtu are not written and the *.vtu file
141 : * contains all information.
142 : *
143 : * ATTENTION: This class uses heavily the mark-function of the grid.
144 : * Do not use any member function while having called begin_mark()
145 : *
146 : * \tparam dim world dimension
147 : */
148 : template <int TDim>
149 : class VTKOutput
150 : {
151 : public:
152 : /// clears the selected output
153 0 : void clear_selection() {m_vSymbFct.clear(); m_vSymbFctNodal.clear(); m_vSymbFctElem.clear();}
154 :
155 : /// clears the selected output
156 0 : void clear_data_selection()
157 : {
158 : m_vScalarNodalData.clear();
159 : m_vScalarElemData.clear();
160 :
161 : m_vVectorNodalData.clear();
162 : m_vVectorElemData.clear();
163 0 : }
164 :
165 : /// schedules that all components of the passed discrete functions will be written to file
166 0 : void select_all(bool bAll) {m_bSelectAll = bAll;}
167 :
168 : /// selects a value of a grid function value to be written
169 : /**
170 : * This function schedules the component(s) passed by symbolic name(s) of an
171 : * approximation space to be written to the vtk file under a specified name.
172 : * The type of the data (nodal/element) will be determined based on the
173 : * trial space of the components (i.e if continuous available or not)
174 : *
175 : * If more than one component is passed, the data will be interpreted as a
176 : * vector and #dim arguments must be passed.
177 : *
178 : * example: fctName = "p"; name = "pressure"
179 : * example: fctNames = "u,v,w"; name = "velocity"
180 : *
181 : * \param[in] fctName symbolic name of component
182 : * \param[in] name name that will appear in the vtk file for the data
183 : */
184 : /// \{
185 : void select(const char* fctName, const char* name);
186 : void select(const std::vector<std::string>& vFct, const char* name);
187 : /// \}
188 :
189 : /// selects a nodal value of a grid function value to be written
190 : /**
191 : * This function schedules the component(s) passed by symbolic name(s) of an
192 : * approximation space to be
193 : * written to the vtk file under a specified name. Note, that for the
194 : * trial space of the component an evaluation of the data at the nodes
195 : * must be available (continuous). If more than one component is passed, the
196 : * data will be interpreted as a vector and #dim arguments must be passed.
197 : *
198 : * example: fctName = "p"; name = "pressure"
199 : * example: fctNames = "u,v,w"; name = "velocity"
200 : *
201 : * \param[in] fctName symbolic name of component
202 : * \param[in] name name that will appear in the vtk file for the data
203 : */
204 : /// \{
205 : void select_nodal(const char* fctName, const char* name);
206 : void select_nodal(const std::vector<std::string>& vFct, const char* name);
207 : /// \}
208 :
209 : /// selects a element value of a grid function value to be written
210 : /**
211 : * This function schedules the component(s) passed by symbolic name(s) of an
212 : * approximation space to be written to the vtk file under a specified name.
213 : * If more than one component is passed, the
214 : * data will be interpreted as a vector and #dim arguments must be passed.
215 : *
216 : * example: fctName = "p"; name = "pressure"
217 : * example: fctNames = "u,v,w"; name = "velocity"
218 : *
219 : * \param[in] fctName symbolic name of component
220 : * \param[in] name name that will appear in the vtk file for the data
221 : */
222 : /// \{
223 : void select_element(const char* fctName, const char* name);
224 : void select_element(const std::vector<std::string>& vFct, const char* name);
225 : /// \}
226 :
227 : /// selects a data value to be written
228 : /**
229 : * This function schedules a user data to be written to the vtk file under
230 : * a specified name.
231 : * The type of the data (nodal/element) will be determined based on the
232 : * trial space of the components (i.e if continuous available or not)
233 : *
234 : * \param[in] spData data to be written
235 : * \param[in] name name that will appear in the vtk file for the data
236 : */
237 : /// \{
238 : void select(SmartPtr<UserData<number, TDim> > spData, const char* name);
239 : void select(SmartPtr<UserData<MathVector<TDim>, TDim> > spData, const char* name);
240 : void select(SmartPtr<UserData<MathMatrix<TDim,TDim>, TDim> > spData, const char* name);
241 : /// \}
242 :
243 : /// selects a nodal data value to be written
244 : /**
245 : * This function schedules a user data to be written to the vtk file under
246 : * a specified name. Note, that for the the data at the nodes must be
247 : * available (continuous).
248 : *
249 : * \param[in] spData data to be written
250 : * \param[in] name name that will appear in the vtk file for the data
251 : */
252 : /// \{
253 : void select_nodal(SmartPtr<UserData<number, TDim> > spData, const char* name);
254 : void select_nodal(SmartPtr<UserData<MathVector<TDim>, TDim> > spData, const char* name);
255 : void select_nodal(SmartPtr<UserData<MathMatrix<TDim,TDim>, TDim> > spData, const char* name);
256 : /// \}
257 :
258 : /// selects a element data value to be written
259 : /**
260 : * This function schedules a user data to be written to the vtk file under
261 : * a specified name.
262 : *
263 : * \param[in] spData data to be written
264 : * \param[in] name name that will appear in the vtk file for the data
265 : */
266 : /// \{
267 : void select_element(SmartPtr<UserData<number, TDim> > spData, const char* name);
268 : void select_element(SmartPtr<UserData<MathVector<TDim>, TDim> > spData, const char* name);
269 : void select_element(SmartPtr<UserData<MathMatrix<TDim,TDim>, TDim> > spData, const char* name);
270 : /// \}
271 :
272 : /**
273 : * This function writes the grid to a *.vtu file. It calls the function
274 : * print_subset for each subset if necessary and produces a grouping
275 : * *.pvd file for paraview. If all subsets have to be written, then only
276 : * a single file is produced containing the whole grid.
277 : *
278 : * \param[in] filename filename for produced files
279 : * \param[in] u grid function
280 : * \param[in] step time step counter (-1 indicates stationary case)
281 : * \param[in] time time point corresponding to timestep
282 : * \param[in] makeConsistent flag if data shall be made consistent before printing
283 : */
284 : template <typename TFunction>
285 : void print(const char* filename, TFunction& u,
286 : int step, number time,
287 : bool makeConsistent);
288 :
289 : /** Calls print with makeConsistent enabled.*/
290 : template <typename TFunction>
291 0 : void print(const char* filename, TFunction& u,
292 : int step, number time)
293 : {
294 0 : return print(filename, u, step, time, true);
295 : }
296 :
297 : /**
298 : * This function simply calls 'print' using step = -1 to indicate the
299 : * stationary case. It is intended to write time independent data.
300 : *
301 : * \param[in] filename filename for produced files
302 : * \param[in] u grid function
303 : * \param[in] makeConsistent flag if data shall be made consistent before printing
304 : */
305 : template <typename TFunction>
306 0 : void print(const char* filename, TFunction& u,
307 : bool makeConsistent)
308 : {
309 0 : print(filename, u, -1, 0.0, makeConsistent);
310 0 : }
311 :
312 : /** Calls print with makeConsistent enabled.*/
313 : template <typename TFunction>
314 0 : void print(const char* filename, TFunction& u)
315 : {
316 0 : return print(filename, u, true);
317 : }
318 :
319 : /// prints the domain to file
320 : void print(const char* filename, Domain<TDim>& domain);
321 :
322 : /**
323 : * This function writes the subset si of the grid (or the whole grid if
324 : * si < 0) to the file "filename.vtu".
325 : *
326 : * If step >= 0 is passed, this indicates that a time step is written and
327 : * to the filename a "*_pXXXX*.vtu" is added, indicating the timestep.
328 : *
329 : * If the computation is done in parallel and the number of Processes is
330 : * greater than one, a "*_pXXXX*.vtu" is added, indicating the process. Then,
331 : * in addition a "filename*.pvtu" file is written, grouping all *vtu files
332 : * from different processes.
333 : *
334 : * If only a part of the grid, i.e. a subset, is written, the to the filename
335 : * a "*_sXXXX*.vtu" is added, to indicate the subset.
336 : *
337 : * \param[in, out] filename base name for output file(s)
338 : * \param[in] u grid function
339 : * \param[in] si Subset (si < 0 indicates whole grid)
340 : * \param[in] step counter for timestep (-1 means stationary)
341 : * \param[in] time time point of timestep
342 : * \param[in] makeConsistent flag if data shall be made consistent before printing
343 : */
344 : template <typename TFunction>
345 : void print_subset(const char* filename, TFunction& u,
346 : int si, int step, number time,
347 : bool makeConsistent);
348 :
349 : /** Calls print_subset with makeConsistent enabled.*/
350 : template <typename TFunction>
351 0 : void print_subset(const char* filename, TFunction& u,
352 : int si, int step = -1, number time = 0.0)
353 : {
354 0 : return print_subset(filename, u, si, step, time, true);
355 : }
356 :
357 : /**
358 : * This function writes the given subsets ssGrp of the grid "filename.vtu".
359 : *
360 : * If step >= 0 is passed, this indicates that a time step is written and
361 : * to the filename a "*_pXXXX*.vtu" is added, indicating the timestep.
362 : *
363 : * If the computation is done in parallel and the number of Processes is
364 : * greater than one, a "*_pXXXX*.vtu" is added, indicating the process. Then,
365 : * in addition a "filename*.pvtu" file is written, grouping all *vtu files
366 : * from different processes.
367 : *
368 : * \param[in, out] filename base name for output file(s)
369 : * \param[in] u grid function
370 : * \param[in] ssGrp subset group
371 : * \param[in] step counter for timestep (-1 means stationary)
372 : * \param[in] time time point of timestep
373 : * \param[in] makeConsistent flag if data shall be made consistent before printing
374 : */
375 : template <typename TFunction>
376 : void print_subsets(const char* filename, TFunction& u,
377 : const SubsetGroup& ssGrp, int step = -1, number time = 0.0,
378 : bool makeConsistent = true);
379 :
380 : /**
381 : * This function writes the given subsets ssGrp of the grid "filename.vtu".
382 : *
383 : * If step >= 0 is passed, this indicates that a time step is written and
384 : * to the filename a "*_pXXXX*.vtu" is added, indicating the timestep.
385 : *
386 : * If the computation is done in parallel and the number of Processes is
387 : * greater than one, a "*_pXXXX*.vtu" is added, indicating the process. Then,
388 : * in addition a "filename*.pvtu" file is written, grouping all *vtu files
389 : * from different processes.
390 : *
391 : * \param[in, out] filename base name for output file(s)
392 : * \param[in] u grid function
393 : * \param[in] ssNames subset names
394 : * \param[in] step counter for timestep (-1 means stationary)
395 : * \param[in] time time point of timestep
396 : * \param[in] makeConsistent flag if data shall be made consistent before printing
397 : */
398 : template <typename TFunction>
399 : void print_subsets(const char* filename, TFunction& u,
400 : const char* ssNames, int step = -1, number time = 0.0,
401 : bool makeConsistent = true);
402 :
403 : /** Calls print_subsets with makeConsistent enabled.*/
404 : template <typename TFunction>
405 0 : void print_subsets(const char* filename, TFunction& u, const char* ssNames,
406 : int step, number time)
407 : {
408 0 : return print_subsets(filename, u, ssNames, step, time, true);
409 : }
410 :
411 : /**
412 : * This function simply calls 'print_subsets' using step = -1 to indicate the
413 : * stationary case. It is intended to write time independent data.
414 : *
415 : * \param[in] filename filename for produced files
416 : * \param[in] u grid function
417 : * \param[in] ssNames subset names
418 : * \param[in] makeConsistent flag if data shall be made consistent before printing
419 : */
420 : template <typename TFunction>
421 0 : void print_subsets(const char* filename, TFunction& u, const char* ssNames,
422 : bool makeConsistent)
423 : {
424 0 : print_subsets(filename, u, ssNames, -1, 0.0, makeConsistent);
425 0 : }
426 :
427 : /** Calls print_subsets with makeConsistent enabled.*/
428 : template <typename TFunction>
429 0 : void print_subsets(const char* filename, TFunction& u, const char* ssNames)
430 : {
431 0 : return print_subsets(filename, u, ssNames, true);
432 : }
433 :
434 : /**
435 : * When a time series has been computed, this function can be used to procduce
436 : * a grouping *.pvd file for paraview visualization.
437 : *
438 : * \param[in] filename filename used in time series
439 : * \param[in] u grid function
440 : */
441 : template <typename TFunction>
442 : void write_time_pvd(const char* filename, TFunction& u);
443 :
444 : /**
445 : * When a time series has been computed, this function can be used to procduce
446 : * a grouping *_processwise.pvd file for paraview visualization.
447 : *
448 : * \param[in] filename filename used in time series
449 : * \param[in] u grid function
450 : */
451 : template <typename TFunction>
452 : void write_time_processwise_pvd(const char* filename, TFunction& u);
453 :
454 : /**
455 : * When a time series has been computed, this function can be used to procduce
456 : * a grouping *.pvd file for paraview visualization.
457 : *
458 : * \param[in] filename filename used in time series
459 : * \param[in] u grid function
460 : * \param[in] si subset index
461 : */
462 : template <typename TFunction>
463 : void write_time_pvd_subset(const char* filename, TFunction& u, int si);
464 : protected:
465 : /**
466 : * This function counts the number of vertices, elements and connections for
467 : * given subsets.
468 : *
469 : * \param[in] u discrete function
470 : * \param[in] ssGrp subset group of the piece
471 : * \param[in] dim dimension of subset
472 : * \param[out] numVert number of vertices
473 : * \param[out] numElem number of elements
474 : * \param[out] numConn number of connections
475 : */
476 : template <typename T>
477 : void
478 : count_piece_sizes(Grid& grid, const T& iterContainer, const SubsetGroup& ssGrp, const int dim,
479 : int& numVert, int& numElem, int& numConn);
480 :
481 : /**
482 : * This function counts the number of vertices, elements and connections that
483 : * are in a subset of the grid (or whole grid, if and only if si < 0).
484 : *
485 : * NOTE: that the number of vertices that are needed to describe the subset
486 : * may be more than the vertices of the subset. This happens, when elements are
487 : * part of the subset, but the vertices of these elements are part of another
488 : * subset. This function counts also those vertices, i.e. it counts the closure
489 : * of all vertices needed to describe the subset.
490 : *
491 : * \param[in] u discrete function
492 : * \param[in] si subset index (si < 0 indicates whole grid)
493 : * \param[out] numVert number of vertices needed to compose subset
494 : * \param[out] numElem number of elements in the subset
495 : * \param[out] numConn number of connections
496 : */
497 : template <typename TElem, typename T>
498 : void
499 : count_sizes(Grid& grid, const T& iterContainer, const int si,
500 : int& numVert, int& numElem, int& numConn);
501 :
502 : template <typename T>
503 : void
504 : write_points_cells_piece(VTKFileWriter& File,
505 : Grid::VertexAttachmentAccessor<Attachment<int> >& aaVrtIndex,
506 : const Grid::VertexAttachmentAccessor<Attachment<MathVector<TDim> > >& aaPos,
507 : Grid& grid, const T& iterContainer, const SubsetGroup& ssGrp,
508 : const int dim, const int numVert, const int numElem, const int numConn);
509 :
510 : template <typename T>
511 : void
512 : write_grid_piece(VTKFileWriter& File,
513 : Grid::VertexAttachmentAccessor<Attachment<int> >& aaVrtIndex,
514 : const Grid::VertexAttachmentAccessor<Attachment<MathVector<TDim> > >& aaPos,
515 : Grid& grid, const T& iterContainer, SubsetGroup& ssGrp, int dim);
516 :
517 : public:
518 : // maybe somebody wants to do this from outside
519 : static void write_empty_grid_piece(VTKFileWriter& File,
520 : bool binary = true);
521 :
522 0 : void set_user_defined_comment(const char* comment) {m_sComment = comment;};
523 :
524 : protected:
525 :
526 : // writes an xml comment to a vtu file
527 : void write_comment(VTKFileWriter& File);
528 :
529 : // writes an xml comment to a pvd file
530 : void write_comment_printf(FILE* File);
531 :
532 : /**
533 : * This function writes the vertices of a piece (the specified subsets) of
534 : * the grid to a vtk file.
535 : *
536 : * \param[in,out] File file to write the points
537 : * \param[in] u discrete function
538 : * \param[in] ssGrp specified subsets
539 : * \param[in] dim dimension of subsets
540 : * \param[in] numVert number of vertices
541 : */
542 : template <typename T>
543 : void
544 : write_points(VTKFileWriter& File,
545 : Grid::VertexAttachmentAccessor<Attachment<int> >& aaVrtIndex,
546 : const Grid::VertexAttachmentAccessor<Attachment<MathVector<TDim> > >& aaPos,
547 : Grid& grid, const T& iterContainer, const SubsetGroup& ssGrp, const int dim,
548 : const int numVert);
549 :
550 : /**
551 : * This method loops all elements of a subset and writes the vertex
552 : * coordinates to a file. All vertices of the element are written reguardless
553 : * if the vertex itself is contained in the subset or another. Written
554 : * vertices are marked and not written again.
555 : *
556 : * \param[in] File file stream
557 : * \param[in] u grid function
558 : * \param[in] si subset index
559 : * \param[in] n counter for vertices
560 : */
561 : template <typename TElem, typename T>
562 : void
563 : write_points_elementwise(VTKFileWriter& File,
564 : Grid::VertexAttachmentAccessor<Attachment<int> >& aaVrtIndex,
565 : const Grid::VertexAttachmentAccessor<Attachment<MathVector<TDim> > >& aaPos,
566 : Grid& grid, const T& iterContainer, const int si, int& n);
567 :
568 : /**
569 : * This function writes the elements that are part of the specified subsets.
570 : *
571 : * \param[in,out] File file to write the points
572 : * \param[in] u discrete function
573 : * \param[in] ssGrp specified subsets
574 : * \param[in] dim dimension of subsets
575 : * \param[out] numElem number of elements
576 : * \param[out] numConn number of connections
577 : */
578 : template <typename T>
579 : void
580 : write_cells(VTKFileWriter& File,
581 : Grid::VertexAttachmentAccessor<Attachment<int> >& aaVrtIndex,
582 : Grid& grid, const T& iterContainer, const SubsetGroup& ssGrp,
583 : const int dim, const int numElem, const int numConn);
584 :
585 :
586 : /**
587 : * This method writes the 'connectivity' for each element of a subset. The
588 : * connectivity are the indices of all vertex the element is formed of.
589 : *
590 : * \param[in] File file stream
591 : * \param[in] u grid function
592 : * \param[in] si subset index
593 : */
594 : template <typename TElem, typename T>
595 : void
596 : write_cell_connectivity(VTKFileWriter& File,
597 : Grid::VertexAttachmentAccessor<Attachment<int> >& aaVrtIndex,
598 : Grid& grid, const T& iterContainer, const int si);
599 :
600 : template <typename T>
601 : void
602 : write_cell_connectivity(VTKFileWriter& File,
603 : Grid::VertexAttachmentAccessor<Attachment<int> >& aaVrtIndex,
604 : Grid& grid, const T& iterContainer, const SubsetGroup& ssGrp,
605 : const int dim, const int numConn);
606 :
607 : /**
608 : * This method writes the 'offset' for each element of a subset.
609 : *
610 : * \param[in] File file stream
611 : * \param[in] u grid function
612 : * \param[in] si subset index
613 : * \param[in] n counter for vertices
614 : */
615 : template <typename TElem, typename T>
616 : void
617 : write_cell_offsets(VTKFileWriter& File, const T& iterContainer, const int si, int& n);
618 :
619 : template <typename T>
620 : void
621 : write_cell_offsets(VTKFileWriter& File, const T& iterContainer, const SubsetGroup& ssGrp,
622 : const int dim, const int numElem);
623 :
624 : /**
625 : * This method writes the 'type' for each element of a subset. The type is
626 : * a vtk-format specified number for a given element type.
627 : *
628 : * \param[in] File file stream
629 : * \param[in] u grid function
630 : * \param[in] si subset index
631 : */
632 : template <typename TElem, typename T>
633 : void
634 : write_cell_types(VTKFileWriter& File, const T& iterContainer, const int si);
635 :
636 : template <typename T>
637 : void
638 : write_cell_types(VTKFileWriter& File, const T& iterContainer, const SubsetGroup& ssGrp,
639 : const int dim, const int numElem);
640 :
641 : /**
642 : * This method writes the 'region' for each element of a subset. The region is
643 : * the ug4 subset index for a given element.
644 : *
645 : * \param[in] File file stream
646 : * \param[in] u grid function
647 : * \param[in] si subset index
648 : * \param[in] sh multigrid subset handler
649 : */
650 : template <typename TElem, typename T>
651 : void
652 : write_cell_subsets(VTKFileWriter& File, const T& iterContainer, const int si, MGSubsetHandler& sh);
653 :
654 : template <typename T>
655 : void
656 : write_cell_subsets(VTKFileWriter& File, const T& iterContainer, const SubsetGroup& ssGrp,
657 : const int dim, const int numElem, MGSubsetHandler& sh);
658 :
659 : /**
660 : * This method writes the subset names.
661 : *
662 : * \param[in] File file stream
663 : * \param[in] sh multigrid subset handler
664 : */
665 : void
666 : write_cell_subset_names(VTKFileWriter& File, MGSubsetHandler& sh);
667 :
668 : /**
669 : * This method writes the proc ranks for each cell.
670 : *
671 : * \param[in] File file stream
672 : * \param[in] u grid function
673 : * \param[in] si subset index
674 : * \param[in] sh multigrid subset handler
675 : */
676 : template <typename TElem, typename T>
677 : void
678 : write_cell_proc_ranks(VTKFileWriter& File, const T& iterContainer, const int si, MGSubsetHandler& sh);
679 :
680 : template <typename T>
681 : void
682 : write_cell_proc_ranks(VTKFileWriter& File, const T& iterContainer, const SubsetGroup& ssGrp,
683 : const int dim, const int numElem, MGSubsetHandler& sh);
684 :
685 : protected:
686 : /**
687 : * This function writes a piece of the grid to the vtk file. First the
688 : * geometric data is written (points, cells, connectivity). Then the data
689 : * associated to the points or cells is written.
690 : *
691 : * \param[in,out] File file to write the points
692 : * \param[in] u discrete function
693 : * \param[in] ssGrp subsets
694 : * \param[in] dim dimension of subset
695 : */
696 : template <typename TFunction>
697 : void
698 : write_grid_solution_piece(VTKFileWriter& File,
699 : Grid::VertexAttachmentAccessor<Attachment<int> >& aaVrtIndex,
700 : Grid& grid, TFunction& u, number time,
701 : const SubsetGroup& ssGrp, const int dim);
702 :
703 : ///////////////////////////////////////////////////////////////////////////
704 : // nodal data
705 :
706 : /**
707 : * This method writes the nodal values of a function to file. If the function
708 : * group consists of more than on component, than a vector of components is
709 : * expected. For each function, access to the (unique) nodal value is
710 : * required.
711 : *
712 : * \param[in] File file stream
713 : * \param[in] u grid function
714 : * \param[in] vFct components to be written
715 : * \param[in] grid Grid
716 : * \param[in] si subset index
717 : */
718 : template <typename TElem, typename TFunction>
719 : void write_nodal_values_elementwise(VTKFileWriter& File, TFunction& u,
720 : const std::vector<size_t>& vFct,
721 : Grid& grid, const int si);
722 :
723 : /**
724 : * This function writes the values of a function as a \<DataArray\> field to
725 : * the vtu file.
726 : *
727 : * \param[in,out] File file to write the points
728 : * \param[in] u discrete function
729 : * \param[in] vFct components to be written
730 : * \param[in] name name to appear for the component
731 : * \param[in] grid Grid
732 : * \param[in] si subset
733 : * \param[in] dim dimension of subset
734 : * \param[in] numVert number of vertices
735 : */
736 : template <typename TFunction>
737 : void write_nodal_values(VTKFileWriter& File, TFunction& u,
738 : const std::vector<size_t>& vFct,
739 : const std::string& name,
740 : Grid& grid, const int si, const int dim,
741 : const int numVert);
742 :
743 : /**
744 : * This function writes the values of a function as a \<DataArray\> field to
745 : * the vtu file.
746 : *
747 : * \param[in,out] File file to write the points
748 : * \param[in] u discrete function
749 : * \param[in] vFct components to be written
750 : * \param[in] name name to appear for the component
751 : * \param[in] grid Grid
752 : * \param[in] ssGrp subsets
753 : * \param[in] dim dimension of subsets
754 : * \param[in] numVert number of vertices
755 : */
756 : template <typename TFunction>
757 : void write_nodal_values(VTKFileWriter& File, TFunction& u,
758 : const std::vector<size_t>& vFct,
759 : const std::string& name,
760 : Grid& grid, const SubsetGroup& ssGrp, const int dim,
761 : const int numVert);
762 :
763 : /// writes the nodal scalar data
764 : /// \{
765 : template <typename TElem, typename TFunction, typename TData>
766 : void write_nodal_data_elementwise(VTKFileWriter& File, TFunction& u,
767 : number time,
768 : SmartPtr<UserData<TData, TDim> > spData,
769 : Grid& grid, const int si);
770 :
771 : template <typename TFunction, typename TData>
772 : void write_nodal_data(VTKFileWriter& File, TFunction& u, number time,
773 : SmartPtr<UserData<TData, TDim> > spData,
774 : const int numCmp,
775 : const std::string& name,
776 : Grid& grid, const SubsetGroup& ssGrp, const int dim,
777 : const int numVert);
778 : /// \}
779 :
780 : template <typename TFunction>
781 : void write_nodal_values_piece(VTKFileWriter& File, TFunction& u, number time,
782 : Grid& grid, const SubsetGroup& ssGrp, const int dim,
783 : const int numVert);
784 :
785 : ///////////////////////////////////////////////////////////////////////////
786 : // cell data
787 :
788 : /**
789 : * This method writes the cell values of a function to file. If the function
790 : * group consists of more than on component, than a vector of components is
791 : * expected.
792 : *
793 : * \param[in] File file stream
794 : * \param[in] u grid function
795 : * \param[in] vFct components to be written
796 : * \param[in] grid Grid
797 : * \param[in] si subset index
798 : */
799 : template <typename TElem, typename TFunction>
800 : void write_cell_values_elementwise(VTKFileWriter& File, TFunction& u,
801 : const std::vector<size_t>& vFct,
802 : Grid& grid, const int si);
803 :
804 : /**
805 : * This function writes the values of a function as a \<DataArray\> field to
806 : * the vtu file.
807 : *
808 : * \param[in,out] File file to write the points
809 : * \param[in] u discrete function
810 : * \param[in] vFct components to be written
811 : * \param[in] name name to appear for the component
812 : * \param[in] grid Grid
813 : * \param[in] ssGrp subset
814 : * \param[in] dim dimension of subsets
815 : * \param[in] numVert number of vertices
816 : */
817 : template <typename TFunction>
818 : void write_cell_values(VTKFileWriter& File, TFunction& u,
819 : const std::vector<size_t>& vFct,
820 : const std::string& name,
821 : Grid& grid, const SubsetGroup& ssGrp, const int dim,
822 : const int numElem);
823 :
824 : /// writes the nodal cell data
825 : /// \{
826 : template <typename TElem, typename TFunction, typename TData>
827 : void write_cell_data_elementwise(VTKFileWriter& File, TFunction& u, number time,
828 : SmartPtr<UserData<TData, TDim> > spData,
829 : Grid& grid, const int si);
830 :
831 : template <typename TFunction, typename TData>
832 : void write_cell_data(VTKFileWriter& File, TFunction& u, number time,
833 : SmartPtr<UserData<TData, TDim> > spData,
834 : const int numCmp,
835 : const std::string& name,
836 : Grid& grid, const SubsetGroup& ssGrp, const int dim,
837 : const int numElem);
838 : /// \}
839 :
840 : template <typename TFunction>
841 : void write_cell_values_piece(VTKFileWriter& File, TFunction& u, number time,
842 : Grid& grid, const SubsetGroup& ssGrp, const int dim,
843 : const int numElem);
844 :
845 : ///////////////////////////////////////////////////////////////////////////
846 : // file names
847 :
848 : /// writes a grouping *.pvtu file
849 : template <typename TFunction>
850 : void write_pvtu(TFunction& u, const std::string& filename,
851 : int si, int step, number time);
852 :
853 : public:
854 : /// writes a grouping *.pvd file, grouping all data from different subsets
855 : static void write_subset_pvd(int numSubset, const std::string& filename,
856 : int step = -1, number time = 0.0);
857 :
858 : /// creates the needed vtu file name
859 : static void vtu_filename(std::string& nameOut, std::string nameIn,
860 : int rank, int si, int maxSi, int step);
861 :
862 : /// create the needed pvtu file name
863 : static void pvtu_filename(std::string& nameOut, std::string nameIn,
864 : int si, int maxSi, int step);
865 :
866 : /// creates the needed pvd file name
867 : static void pvd_filename(std::string& nameOut, std::string nameIn);
868 :
869 : /// creates the needed pvd file name to group the time steps
870 : static void pvd_time_filename(std::string& nameOut, std::string nameIn,
871 : int step);
872 :
873 : public:
874 : /// default constructor
875 0 : VTKOutput() : m_bSelectAll(true), m_bBinary(true), m_bWriteGrid(true), m_bWriteSubsetIndices(false), m_bWriteProcRanks(false) {} //TODO: maybe true?
876 :
877 : /// should values be printed in binary (base64 encoded way ) or plain ascii
878 0 : void set_binary(bool b) {m_bBinary = b;};
879 :
880 0 : void set_write_grid(bool b) {m_bWriteGrid = b;};
881 :
882 0 : void set_write_subset_indices(bool b) {m_bWriteSubsetIndices = b;};
883 :
884 0 : void set_write_proc_ranks(bool b) {m_bWriteProcRanks = b;};
885 :
886 : protected:
887 : /// returns true if name for vtk-component is already used
888 : bool vtk_name_used(const char* name) const;
889 :
890 : /// writes data to stream
891 : /**
892 : * The purpose of the function is to convert a double data to binary float
893 : * (the Float32-format as used in vtk) or ASCII representation. Note that
894 : * the ASCII output should contain only numbers in the range of the normalized
895 : * Float32 representation, i.e. with the minimum absolute value
896 : * 1.1755e-38 = 2^-126. As 'float'-variables can potentially store denormalized
897 : * numbers (and this is really the case in the practice) with the minimum
898 : * absolute value 1.4013e-45 = 2^-149 that are not correctly read by some
899 : * visualization tools (in some operating systems), these function replace
900 : * those values (lying near to 0) with 0.
901 : */
902 : /// \{
903 : inline void write_item_to_file(VTKFileWriter& File, float data);
904 0 : inline void write_item_to_file(VTKFileWriter& File, double data) {write_item_to_file(File, (float) data);};
905 : inline void write_item_to_file(VTKFileWriter& File, const ug::MathVector<1>& data);
906 : inline void write_item_to_file(VTKFileWriter& File, const ug::MathVector<2>& data);
907 : inline void write_item_to_file(VTKFileWriter& File, const ug::MathVector<3>& data);
908 : inline void write_item_to_file(VTKFileWriter& File, const ug::MathMatrix<1,1>& data);
909 : inline void write_item_to_file(VTKFileWriter& File, const ug::MathMatrix<2,2>& data);
910 : inline void write_item_to_file(VTKFileWriter& File, const ug::MathMatrix<3,3>& data);
911 : /// \}
912 :
913 : /// prints ascii representation of a float in the Float32 format (a protection against the denormalized floats)
914 : inline VTKFileWriter& write_asc_float(VTKFileWriter& File, float data);
915 :
916 : protected:
917 : /// scheduled components to be printed
918 : bool m_bSelectAll;
919 : /// print values in binary (base64 encoded way) or plain ascii
920 : bool m_bBinary;
921 : std::map<std::string, std::vector<std::string> > m_vSymbFct;
922 : std::map<std::string, std::vector<std::string> > m_vSymbFctNodal;
923 : std::map<std::string, std::vector<std::string> > m_vSymbFctElem;
924 : typedef typename std::map<std::string,
925 : std::vector<std::string> >::iterator ComponentsIterator;
926 :
927 : /// scheduled scalar data to be printed
928 : std::map<std::string, SmartPtr<UserData<number, TDim> > > m_vScalarNodalData;
929 : std::map<std::string, SmartPtr<UserData<number, TDim> > > m_vScalarElemData;
930 : typedef typename std::map<std::string,
931 : SmartPtr<UserData<number, TDim> > >::iterator ScalarDataIterator;
932 :
933 : /// scheduled vector data to be printed
934 : std::map<std::string, SmartPtr<UserData<MathVector<TDim>, TDim> > > m_vVectorNodalData;
935 : std::map<std::string, SmartPtr<UserData<MathVector<TDim>, TDim> > > m_vVectorElemData;
936 : typedef typename std::map<std::string,
937 : SmartPtr<UserData<MathVector<TDim>, TDim> > >::iterator VectorDataIterator;
938 :
939 : /// scheduled vector data to be printed
940 : std::map<std::string, SmartPtr<UserData<MathMatrix<TDim,TDim>, TDim> > > m_vMatrixNodalData;
941 : std::map<std::string, SmartPtr<UserData<MathMatrix<TDim,TDim>, TDim> > > m_vMatrixElemData;
942 : typedef typename std::map<std::string,
943 : SmartPtr<UserData<MathMatrix<TDim,TDim>, TDim> > >::iterator MatrixDataIterator;
944 :
945 : /// map storing the time points
946 : std::map<std::string, std::vector<number> > m_mTimestep;
947 :
948 : bool m_bWriteGrid;
949 : std::string m_sComment;
950 :
951 : bool m_bWriteSubsetIndices;
952 : bool m_bWriteProcRanks;
953 : };
954 :
955 : } // namespace ug
956 :
957 : #include "vtkoutput_impl.h"
958 :
959 : #endif /* __H__UG__LIB_DISC__IO__VTKOUTPUT__ */
|