ug4
vtkoutput.h
Go to the documentation of this file.
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
47 #include "lib_disc/domain.h"
49 
50 namespace ug{
51 // todo this should avoid refactoring all the signatures of VTKOutput, remove it later
53 
54 template <typename T>
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  static typename traits<TElem>::const_iterator begin(const T& provider, int si)
66  {
67  if(si < 0) return provider.template begin<TElem>();
68  else return provider.template begin<TElem>(si);
69  }
70 
71  template <typename TElem>
72  static typename traits<TElem>::const_iterator end(const T& provider, int si)
73  {
74  if(si < 0) return provider.template end<TElem>();
75  else return provider.template end<TElem>(si);
76  }
77 };
78 
79 template <>
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  static typename traits<TElem>::const_iterator begin(const T& provider, int si)
95  {
96  const int lev = provider.num_levels() - 1;
97  if(si < 0) return provider.multi_grid()->template begin<TElem>(lev);
98  else return provider.template begin<TElem>(si, lev);
99  }
100 
101  template <typename TElem>
102  static typename traits<TElem>::const_iterator end(const T& provider, int si)
103  {
104  const int lev = provider.num_levels() - 1;
105  if(si < 0) return provider.multi_grid()->template end<TElem>(lev);
106  else return provider.template end<TElem>(si, lev);
107  }
108 };
109 
110 
112 
148 template <int TDim>
150 {
151  public:
153  void clear_selection() {m_vSymbFct.clear(); m_vSymbFctNodal.clear(); m_vSymbFctElem.clear();}
154 
157  {
158  m_vScalarNodalData.clear();
159  m_vScalarElemData.clear();
160 
161  m_vVectorNodalData.clear();
162  m_vVectorElemData.clear();
163  }
164 
166  void select_all(bool bAll) {m_bSelectAll = bAll;}
167 
169 
185  void select(const char* fctName, const char* name);
186  void select(const std::vector<std::string>& vFct, const char* name);
188 
190 
205  void select_nodal(const char* fctName, const char* name);
206  void select_nodal(const std::vector<std::string>& vFct, const char* name);
208 
210 
223  void select_element(const char* fctName, const char* name);
224  void select_element(const std::vector<std::string>& vFct, const char* name);
226 
228 
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);
242 
244 
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);
257 
259 
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);
271 
284  template <typename TFunction>
285  void print(const char* filename, TFunction& u,
286  int step, number time,
287  bool makeConsistent);
288 
290  template <typename TFunction>
291  void print(const char* filename, TFunction& u,
292  int step, number time)
293  {
294  return print(filename, u, step, time, true);
295  }
296 
305  template <typename TFunction>
306  void print(const char* filename, TFunction& u,
307  bool makeConsistent)
308  {
309  print(filename, u, -1, 0.0, makeConsistent);
310  }
311 
313  template <typename TFunction>
314  void print(const char* filename, TFunction& u)
315  {
316  return print(filename, u, true);
317  }
318 
320  void print(const char* filename, Domain<TDim>& domain);
321 
344  template <typename TFunction>
345  void print_subset(const char* filename, TFunction& u,
346  int si, int step, number time,
347  bool makeConsistent);
348 
350  template <typename TFunction>
351  void print_subset(const char* filename, TFunction& u,
352  int si, int step = -1, number time = 0.0)
353  {
354  return print_subset(filename, u, si, step, time, true);
355  }
356 
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 
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 
404  template <typename TFunction>
405  void print_subsets(const char* filename, TFunction& u, const char* ssNames,
406  int step, number time)
407  {
408  return print_subsets(filename, u, ssNames, step, time, true);
409  }
410 
420  template <typename TFunction>
421  void print_subsets(const char* filename, TFunction& u, const char* ssNames,
422  bool makeConsistent)
423  {
424  print_subsets(filename, u, ssNames, -1, 0.0, makeConsistent);
425  }
426 
428  template <typename TFunction>
429  void print_subsets(const char* filename, TFunction& u, const char* ssNames)
430  {
431  return print_subsets(filename, u, ssNames, true);
432  }
433 
441  template <typename TFunction>
442  void write_time_pvd(const char* filename, TFunction& u);
443 
451  template <typename TFunction>
452  void write_time_processwise_pvd(const char* filename, TFunction& u);
453 
462  template <typename TFunction>
463  void write_time_pvd_subset(const char* filename, TFunction& u, int si);
464 protected:
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 
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
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
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  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 
542  template <typename T>
543  void
547  Grid& grid, const T& iterContainer, const SubsetGroup& ssGrp, const int dim,
548  const int numVert);
549 
561  template <typename TElem, typename T>
562  void
566  Grid& grid, const T& iterContainer, const int si, int& n);
567 
578  template <typename T>
579  void
582  Grid& grid, const T& iterContainer, const SubsetGroup& ssGrp,
583  const int dim, const int numElem, const int numConn);
584 
585 
594  template <typename TElem, typename T>
595  void
598  Grid& grid, const T& iterContainer, const int si);
599 
600  template <typename T>
601  void
604  Grid& grid, const T& iterContainer, const SubsetGroup& ssGrp,
605  const int dim, const int numConn);
606 
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 
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 
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 
665  void
667 
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:
696  template <typename TFunction>
697  void
700  Grid& grid, TFunction& u, number time,
701  const SubsetGroup& ssGrp, const int dim);
702 
704  // nodal data
705 
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 
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 
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 
765  template <typename TElem, typename TFunction, typename TData>
766  void write_nodal_data_elementwise(VTKFileWriter& File, TFunction& u,
767  number time,
769  Grid& grid, const int si);
770 
771  template <typename TFunction, typename TData>
772  void write_nodal_data(VTKFileWriter& File, TFunction& u, number time,
774  const int numCmp,
775  const std::string& name,
776  Grid& grid, const SubsetGroup& ssGrp, const int dim,
777  const int numVert);
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 
786  // cell data
787 
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 
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 
826  template <typename TElem, typename TFunction, typename TData>
827  void write_cell_data_elementwise(VTKFileWriter& File, TFunction& u, number time,
829  Grid& grid, const int si);
830 
831  template <typename TFunction, typename TData>
832  void write_cell_data(VTKFileWriter& File, TFunction& u, number time,
834  const int numCmp,
835  const std::string& name,
836  Grid& grid, const SubsetGroup& ssGrp, const int dim,
837  const int numElem);
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 
846  // file names
847 
849  template <typename TFunction>
850  void write_pvtu(TFunction& u, const std::string& filename,
851  int si, int step, number time);
852 
853  public:
855  static void write_subset_pvd(int numSubset, const std::string& filename,
856  int step = -1, number time = 0.0);
857 
859  static void vtu_filename(std::string& nameOut, std::string nameIn,
860  int rank, int si, int maxSi, int step);
861 
863  static void pvtu_filename(std::string& nameOut, std::string nameIn,
864  int si, int maxSi, int step);
865 
867  static void pvd_filename(std::string& nameOut, std::string nameIn);
868 
870  static void pvd_time_filename(std::string& nameOut, std::string nameIn,
871  int step);
872 
873  public:
875  VTKOutput() : m_bSelectAll(true), m_bBinary(true), m_bWriteGrid(true), m_bWriteSubsetIndices(false), m_bWriteProcRanks(false) {} //TODO: maybe true?
876 
878  void set_binary(bool b) {m_bBinary = b;};
879 
880  void set_write_grid(bool b) {m_bWriteGrid = b;};
881 
883 
885 
886  protected:
888  bool vtk_name_used(const char* name) const;
889 
891 
903  inline void write_item_to_file(VTKFileWriter& File, float data);
904  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);
912 
914  inline VTKFileWriter& write_asc_float(VTKFileWriter& File, float data);
915 
916  protected:
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 
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,
932 
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,
938 
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,
944 
946  std::map<std::string, std::vector<number> > m_mTimestep;
947 
949  std::string m_sComment;
950 
953 };
954 
955 } // namespace ug
956 
957 #include "vtkoutput_impl.h"
958 
959 #endif /* __H__UG__LIB_DISC__IO__VTKOUTPUT__ */
location name
Definition: checkpoint_util.lua:128
Definition: smart_pointer.h:108
File writer allowing selective base64 encoding of arbitrary data.
Definition: base64_file_writer.h:75
Definition: domain.h:287
Manages the elements of a grid and their interconnection.
Definition: grid.h:132
A class for fixed size, dense matrices.
Definition: math_matrix.h:52
a mathematical Vector with N entries.
Definition: math_vector.h:97
Handles subsets on a per level basis.
Definition: subset_handler_multi_grid.h:60
MultiGrid * multi_grid()
Definition: subset_handler_multi_grid.h:73
uint num_levels() const
returns the number of levels
Definition: subset_handler_multi_grid.h:80
Group of subsets.
Definition: subset_group.h:51
Type based UserData.
Definition: user_data.h:117
output writer to the VTK file format
Definition: vtkoutput.h:150
void print(const char *filename, TFunction &u, bool makeConsistent)
Definition: vtkoutput.h:306
void count_sizes(Grid &grid, const T &iterContainer, const int si, int &numVert, int &numElem, int &numConn)
Definition: vtkoutput_impl.h:655
void select(const char *fctName, const char *name)
selects a value of a grid function value to be written
Definition: vtkoutput.cpp:440
void write_nodal_data_elementwise(VTKFileWriter &File, TFunction &u, number time, SmartPtr< UserData< TData, TDim > > spData, Grid &grid, const int si)
Definition: vtkoutput_impl.h:1441
std::map< std::string, SmartPtr< UserData< MathVector< TDim >, TDim > > >::iterator VectorDataIterator
Definition: vtkoutput.h:937
void write_time_processwise_pvd(const char *filename, TFunction &u)
Definition: vtkoutput_impl.h:2536
void print_subset(const char *filename, TFunction &u, int si, int step=-1, number time=0.0)
Definition: vtkoutput.h:351
void write_cell_values(VTKFileWriter &File, TFunction &u, const std::vector< size_t > &vFct, const std::string &name, Grid &grid, const SubsetGroup &ssGrp, const int dim, const int numElem)
Definition: vtkoutput_impl.h:2069
std::map< std::string, std::vector< std::string > >::iterator ComponentsIterator
Definition: vtkoutput.h:925
void write_item_to_file(VTKFileWriter &File, float data)
writes data to stream
Definition: vtkoutput_impl.h:87
void write_pvtu(TFunction &u, const std::string &filename, int si, int step, number time)
writes a grouping *.pvtu file
Definition: vtkoutput_impl.h:2230
void write_nodal_values_elementwise(VTKFileWriter &File, TFunction &u, const std::vector< size_t > &vFct, Grid &grid, const int si)
Definition: vtkoutput_impl.h:1593
void set_write_subset_indices(bool b)
Definition: vtkoutput.h:882
void write_nodal_values_piece(VTKFileWriter &File, TFunction &u, number time, Grid &grid, const SubsetGroup &ssGrp, const int dim, const int numVert)
Definition: vtkoutput_impl.h:1710
void write_cell_offsets(VTKFileWriter &File, const T &iterContainer, const int si, int &n)
Definition: vtkoutput_impl.h:1059
void write_time_pvd(const char *filename, TFunction &u)
Definition: vtkoutput_impl.h:2441
bool m_bSelectAll
scheduled components to be printed
Definition: vtkoutput.h:918
std::map< std::string, SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > > > m_vMatrixElemData
Definition: vtkoutput.h:941
std::map< std::string, SmartPtr< UserData< number, TDim > > > m_vScalarElemData
Definition: vtkoutput.h:929
void write_cell_subsets(VTKFileWriter &File, const T &iterContainer, const int si, MGSubsetHandler &sh)
Definition: vtkoutput_impl.h:1277
std::map< std::string, std::vector< std::string > > m_vSymbFctElem
Definition: vtkoutput.h:923
bool vtk_name_used(const char *name) const
returns true if name for vtk-component is already used
Definition: vtkoutput.cpp:665
void clear_selection()
clears the selected output
Definition: vtkoutput.h:153
void write_grid_piece(VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, const Grid::VertexAttachmentAccessor< Attachment< MathVector< TDim > > > &aaPos, Grid &grid, const T &iterContainer, SubsetGroup &ssGrp, int dim)
Definition: vtkoutput_impl.h:563
void write_cell_values_piece(VTKFileWriter &File, TFunction &u, number time, Grid &grid, const SubsetGroup &ssGrp, const int dim, const int numElem)
Definition: vtkoutput_impl.h:2110
void write_item_to_file(VTKFileWriter &File, double data)
writes data to stream
Definition: vtkoutput.h:904
void write_cell_values_elementwise(VTKFileWriter &File, TFunction &u, const std::vector< size_t > &vFct, Grid &grid, const int si)
Definition: vtkoutput_impl.h:1973
void select_all(bool bAll)
schedules that all components of the passed discrete functions will be written to file
Definition: vtkoutput.h:166
void print_subsets(const char *filename, TFunction &u, const SubsetGroup &ssGrp, int step=-1, number time=0.0, bool makeConsistent=true)
Definition: vtkoutput_impl.h:382
void write_nodal_data(VTKFileWriter &File, TFunction &u, number time, SmartPtr< UserData< TData, TDim > > spData, const int numCmp, const std::string &name, Grid &grid, const SubsetGroup &ssGrp, const int dim, const int numVert)
Definition: vtkoutput_impl.h:1536
void count_piece_sizes(Grid &grid, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, int &numVert, int &numElem, int &numConn)
Definition: vtkoutput_impl.h:713
void set_write_proc_ranks(bool b)
Definition: vtkoutput.h:884
void write_comment_printf(FILE *File)
Definition: vtkoutput.cpp:193
std::map< std::string, std::vector< std::string > > m_vSymbFctNodal
Definition: vtkoutput.h:922
void print_subsets(const char *filename, TFunction &u, const char *ssNames, bool makeConsistent)
Definition: vtkoutput.h:421
std::map< std::string, SmartPtr< UserData< number, TDim > > > m_vScalarNodalData
scheduled scalar data to be printed
Definition: vtkoutput.h:928
std::map< std::string, SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > > > m_vMatrixNodalData
scheduled vector data to be printed
Definition: vtkoutput.h:940
std::map< std::string, SmartPtr< UserData< MathVector< TDim >, TDim > > > m_vVectorNodalData
scheduled vector data to be printed
Definition: vtkoutput.h:934
static void write_subset_pvd(int numSubset, const std::string &filename, int step=-1, number time=0.0)
writes a grouping *.pvd file, grouping all data from different subsets
Definition: vtkoutput.cpp:307
void select_element(const char *fctName, const char *name)
selects a element value of a grid function value to be written
Definition: vtkoutput.cpp:607
void print(const char *filename, TFunction &u)
Definition: vtkoutput.h:314
VTKOutput()
default constructor
Definition: vtkoutput.h:875
void write_cell_data_elementwise(VTKFileWriter &File, TFunction &u, number time, SmartPtr< UserData< TData, TDim > > spData, Grid &grid, const int si)
Definition: vtkoutput_impl.h:1831
void clear_data_selection()
clears the selected output
Definition: vtkoutput.h:156
void set_binary(bool b)
should values be printed in binary (base64 encoded way ) or plain ascii
Definition: vtkoutput.h:878
void set_write_grid(bool b)
Definition: vtkoutput.h:880
static void pvd_filename(std::string &nameOut, std::string nameIn)
creates the needed pvd file name
Definition: vtkoutput.cpp:280
void write_time_pvd_subset(const char *filename, TFunction &u, int si)
Definition: vtkoutput_impl.h:2611
static void pvd_time_filename(std::string &nameOut, std::string nameIn, int step)
creates the needed pvd file name to group the time steps
Definition: vtkoutput.cpp:292
void write_nodal_values(VTKFileWriter &File, TFunction &u, const std::vector< size_t > &vFct, const std::string &name, Grid &grid, const int si, const int dim, const int numVert)
static void write_empty_grid_piece(VTKFileWriter &File, bool binary=true)
Definition: vtkoutput.cpp:136
std::map< std::string, std::vector< std::string > > m_vSymbFct
Definition: vtkoutput.h:921
std::map< std::string, std::vector< number > > m_mTimestep
map storing the time points
Definition: vtkoutput.h:946
static void pvtu_filename(std::string &nameOut, std::string nameIn, int si, int maxSi, int step)
create the needed pvtu file name
Definition: vtkoutput.cpp:259
void set_user_defined_comment(const char *comment)
Definition: vtkoutput.h:522
void write_cell_connectivity(VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, Grid &grid, const T &iterContainer, const int si)
Definition: vtkoutput_impl.h:920
void print_subset(const char *filename, TFunction &u, int si, int step, number time, bool makeConsistent)
Definition: vtkoutput_impl.h:247
void print_subsets(const char *filename, TFunction &u, const char *ssNames)
Definition: vtkoutput.h:429
bool m_bWriteProcRanks
Definition: vtkoutput.h:952
void select_nodal(const char *fctName, const char *name)
selects a nodal value of a grid function value to be written
Definition: vtkoutput.cpp:515
std::map< std::string, SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > > >::iterator MatrixDataIterator
Definition: vtkoutput.h:943
bool m_bWriteGrid
Definition: vtkoutput.h:948
void print(const char *filename, TFunction &u, int step, number time, bool makeConsistent)
Definition: vtkoutput_impl.h:184
bool m_bBinary
print values in binary (base64 encoded way) or plain ascii
Definition: vtkoutput.h:920
void print(const char *filename, TFunction &u, int step, number time)
Definition: vtkoutput.h:291
void write_points_elementwise(VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, const Grid::VertexAttachmentAccessor< Attachment< MathVector< TDim > > > &aaPos, Grid &grid, const T &iterContainer, const int si, int &n)
Definition: vtkoutput_impl.h:760
void print_subsets(const char *filename, TFunction &u, const char *ssNames, int step, number time)
Definition: vtkoutput.h:405
void write_points(VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, const Grid::VertexAttachmentAccessor< Attachment< MathVector< TDim > > > &aaPos, Grid &grid, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, const int numVert)
Definition: vtkoutput_impl.h:811
static void vtu_filename(std::string &nameOut, std::string nameIn, int rank, int si, int maxSi, int step)
creates the needed vtu file name
Definition: vtkoutput.cpp:232
std::map< std::string, SmartPtr< UserData< MathVector< TDim >, TDim > > > m_vVectorElemData
Definition: vtkoutput.h:935
void write_cells(VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, Grid &grid, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, const int numElem, const int numConn)
Definition: vtkoutput_impl.h:888
void write_cell_subset_names(VTKFileWriter &File, MGSubsetHandler &sh)
Definition: vtkoutput_impl.h:872
std::string m_sComment
Definition: vtkoutput.h:949
void write_cell_proc_ranks(VTKFileWriter &File, const T &iterContainer, const int si, MGSubsetHandler &sh)
Definition: vtkoutput_impl.h:1359
bool m_bWriteSubsetIndices
Definition: vtkoutput.h:951
void write_cell_data(VTKFileWriter &File, TFunction &u, number time, SmartPtr< UserData< TData, TDim > > spData, const int numCmp, const std::string &name, Grid &grid, const SubsetGroup &ssGrp, const int dim, const int numElem)
Definition: vtkoutput_impl.h:1927
void write_points_cells_piece(VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, const Grid::VertexAttachmentAccessor< Attachment< MathVector< TDim > > > &aaPos, Grid &grid, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, const int numVert, const int numElem, const int numConn)
Definition: vtkoutput_impl.h:538
void write_comment(VTKFileWriter &File)
Definition: vtkoutput.cpp:179
std::map< std::string, SmartPtr< UserData< number, TDim > > >::iterator ScalarDataIterator
Definition: vtkoutput.h:931
VTKFileWriter & write_asc_float(VTKFileWriter &File, float data)
prints ascii representation of a float in the Float32 format (a protection against the denormalized f...
Definition: vtkoutput_impl.h:69
void write_grid_solution_piece(VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, Grid &grid, TFunction &u, number time, const SubsetGroup &ssGrp, const int dim)
Definition: vtkoutput_impl.h:600
void write_cell_types(VTKFileWriter &File, const T &iterContainer, const int si)
Definition: vtkoutput_impl.h:1163
static const int dim
double number
Definition: types.h:124
the ug namespace
Base64FileWriter VTKFileWriter
Definition: vtkoutput.h:52
Definition: vtkoutput.h:59
T::template traits< TElem >::const_iterator const_iterator
Definition: vtkoutput.h:61
T::template traits< TElem >::iterator iterator
Definition: vtkoutput.h:60
T::template traits< TElem >::const_iterator const_iterator
Definition: vtkoutput.h:90
T::template traits< TElem >::iterator iterator
Definition: vtkoutput.h:89
MGSubsetHandler T
Definition: vtkoutput.h:83
static traits< TElem >::const_iterator begin(const T &provider, int si)
Definition: vtkoutput.h:94
static traits< TElem >::const_iterator end(const T &provider, int si)
Definition: vtkoutput.h:102
Definition: vtkoutput.h:56
static traits< TElem >::const_iterator begin(const T &provider, int si)
Definition: vtkoutput.h:65
static traits< TElem >::const_iterator end(const T &provider, int si)
Definition: vtkoutput.h:72