ug4
ug::VTKOutput< TDim > Class Template Reference

output writer to the VTK file format More...

#include <vtkoutput.h>

Public Member Functions

void clear_data_selection ()
 clears the selected output More...
 
void clear_selection ()
 clears the selected output More...
 
void print (const char *filename, Domain< TDim > &domain)
 prints the domain to file More...
 
template<typename TFunction >
void print (const char *filename, TFunction &u)
 
template<typename TFunction >
void print (const char *filename, TFunction &u, bool makeConsistent)
 
template<typename TFunction >
void print (const char *filename, TFunction &u, int step, number time)
 
template<typename TFunction >
void print (const char *filename, TFunction &u, int step, number time, bool makeConsistent)
 
template<typename TFunction >
void print_subset (const char *filename, TFunction &u, int si, int step, number time, bool makeConsistent)
 
template<typename TFunction >
void print_subset (const char *filename, TFunction &u, int si, int step=-1, number time=0.0)
 
template<typename TFunction >
void print_subsets (const char *filename, TFunction &u, const char *ssNames)
 
template<typename TFunction >
void print_subsets (const char *filename, TFunction &u, const char *ssNames, bool makeConsistent)
 
template<typename TFunction >
void print_subsets (const char *filename, TFunction &u, const char *ssNames, int step, number time)
 
template<typename TFunction >
void print_subsets (const char *filename, TFunction &u, const char *ssNames, int step=-1, number time=0.0, bool makeConsistent=true)
 
template<typename TFunction >
void print_subsets (const char *filename, TFunction &u, const SubsetGroup &ssGrp, int step=-1, number time=0.0, bool makeConsistent=true)
 
void select_all (bool bAll)
 schedules that all components of the passed discrete functions will be written to file More...
 
void set_binary (bool b)
 should values be printed in binary (base64 encoded way ) or plain ascii More...
 
void set_user_defined_comment (const char *comment)
 
void set_write_grid (bool b)
 
void set_write_proc_ranks (bool b)
 
void set_write_subset_indices (bool b)
 
 VTKOutput ()
 default constructor More...
 
template<typename TFunction >
void write_time_processwise_pvd (const char *filename, TFunction &u)
 
template<typename TFunction >
void write_time_pvd (const char *filename, TFunction &u)
 
template<typename TFunction >
void write_time_pvd_subset (const char *filename, TFunction &u, int si)
 
void select (const char *fctName, const char *name)
 selects a value of a grid function value to be written More...
 
void select (const std::vector< std::string > &vFct, const char *name)
 selects a value of a grid function value to be written More...
 
void select_nodal (const char *fctName, const char *name)
 selects a nodal value of a grid function value to be written More...
 
void select_nodal (const std::vector< std::string > &vFct, const char *name)
 selects a nodal value of a grid function value to be written More...
 
void select_element (const char *fctName, const char *name)
 selects a element value of a grid function value to be written More...
 
void select_element (const std::vector< std::string > &vFct, const char *name)
 selects a element value of a grid function value to be written More...
 
void select (SmartPtr< UserData< number, TDim > > spData, const char *name)
 selects a data value to be written More...
 
void select (SmartPtr< UserData< MathVector< TDim >, TDim > > spData, const char *name)
 selects a data value to be written More...
 
void select (SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > > spData, const char *name)
 selects a data value to be written More...
 
void select_nodal (SmartPtr< UserData< number, TDim > > spData, const char *name)
 selects a nodal data value to be written More...
 
void select_nodal (SmartPtr< UserData< MathVector< TDim >, TDim > > spData, const char *name)
 selects a nodal data value to be written More...
 
void select_nodal (SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > > spData, const char *name)
 selects a nodal data value to be written More...
 
void select_element (SmartPtr< UserData< number, TDim > > spData, const char *name)
 selects a element data value to be written More...
 
void select_element (SmartPtr< UserData< MathVector< TDim >, TDim > > spData, const char *name)
 selects a element data value to be written More...
 
void select_element (SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > > spData, const char *name)
 selects a element data value to be written More...
 

Static Public Member Functions

static void pvd_filename (std::string &nameOut, std::string nameIn)
 creates the needed pvd file name More...
 
static void pvd_time_filename (std::string &nameOut, std::string nameIn, int step)
 creates the needed pvd file name to group the time steps More...
 
static void pvtu_filename (std::string &nameOut, std::string nameIn, int si, int maxSi, int step)
 create the needed pvtu file name More...
 
static void vtu_filename (std::string &nameOut, std::string nameIn, int rank, int si, int maxSi, int step)
 creates the needed vtu file name More...
 
static void write_empty_grid_piece (VTKFileWriter &File, bool binary=true)
 
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 More...
 

Protected Types

typedef std::map< std::string, std::vector< std::string > >::iterator ComponentsIterator
 
typedef std::map< std::string, SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > > >::iterator MatrixDataIterator
 
typedef std::map< std::string, SmartPtr< UserData< number, TDim > > >::iterator ScalarDataIterator
 
typedef std::map< std::string, SmartPtr< UserData< MathVector< TDim >, TDim > > >::iterator VectorDataIterator
 

Protected Member Functions

template<typename T >
void count_piece_sizes (Grid &grid, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, int &numVert, int &numElem, int &numConn)
 
template<typename TElem , typename T >
void count_sizes (Grid &grid, const T &iterContainer, const int si, int &numVert, int &numElem, int &numConn)
 
bool vtk_name_used (const char *name) const
 returns true if name for vtk-component is already used More...
 
VTKFileWriterwrite_asc_float (VTKFileWriter &File, float data)
 prints ascii representation of a float in the Float32 format (a protection against the denormalized floats) More...
 
template<typename TElem , typename T >
void write_cell_connectivity (VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, Grid &grid, const T &iterContainer, const int si)
 
template<typename T >
void write_cell_connectivity (VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, Grid &grid, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, const int numConn)
 
template<typename TElem , typename T >
void write_cell_offsets (VTKFileWriter &File, const T &iterContainer, const int si, int &n)
 
template<typename T >
void write_cell_offsets (VTKFileWriter &File, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, const int numElem)
 
template<typename TElem , typename T >
void write_cell_proc_ranks (VTKFileWriter &File, const T &iterContainer, const int si, MGSubsetHandler &sh)
 
template<typename T >
void write_cell_proc_ranks (VTKFileWriter &File, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, const int numElem, MGSubsetHandler &sh)
 
void write_cell_subset_names (VTKFileWriter &File, MGSubsetHandler &sh)
 
template<typename TElem , typename T >
void write_cell_subsets (VTKFileWriter &File, const T &iterContainer, const int si, MGSubsetHandler &sh)
 
template<typename T >
void write_cell_subsets (VTKFileWriter &File, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, const int numElem, MGSubsetHandler &sh)
 
template<typename TElem , typename T >
void write_cell_types (VTKFileWriter &File, const T &iterContainer, const int si)
 
template<typename T >
void write_cell_types (VTKFileWriter &File, const T &iterContainer, const SubsetGroup &ssGrp, const int dim, const int numElem)
 
template<typename TFunction >
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)
 
template<typename TElem , typename TFunction >
void write_cell_values_elementwise (VTKFileWriter &File, TFunction &u, const std::vector< size_t > &vFct, Grid &grid, const int si)
 
template<typename TFunction >
void write_cell_values_piece (VTKFileWriter &File, TFunction &u, number time, Grid &grid, const SubsetGroup &ssGrp, const int dim, const int numElem)
 
template<typename T >
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)
 
void write_comment (VTKFileWriter &File)
 
void write_comment_printf (FILE *File)
 
template<typename T >
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)
 
template<typename TFunction >
void write_grid_solution_piece (VTKFileWriter &File, Grid::VertexAttachmentAccessor< Attachment< int > > &aaVrtIndex, Grid &grid, TFunction &u, number time, const SubsetGroup &ssGrp, const int dim)
 
template<typename TFunction >
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)
 
template<typename TFunction >
void write_nodal_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 numVert)
 
template<typename TElem , typename TFunction >
void write_nodal_values_elementwise (VTKFileWriter &File, TFunction &u, const std::vector< size_t > &vFct, Grid &grid, const int si)
 
template<typename TFunction >
void write_nodal_values_piece (VTKFileWriter &File, TFunction &u, number time, Grid &grid, const SubsetGroup &ssGrp, const int dim, const int numVert)
 
template<typename T >
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)
 
template<typename T >
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)
 
template<typename TElem , typename T >
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)
 
template<typename TFunction >
void write_pvtu (TFunction &u, const std::string &filename, int si, int step, number time)
 writes a grouping *.pvtu file More...
 
template<typename TElem , typename TFunction , typename TData >
void write_nodal_data_elementwise (VTKFileWriter &File, TFunction &u, number time, SmartPtr< UserData< TData, TDim > > spData, Grid &grid, const int si)
 
template<typename TFunction , typename TData >
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)
 
template<typename TElem , typename TFunction , typename TData >
void write_cell_data_elementwise (VTKFileWriter &File, TFunction &u, number time, SmartPtr< UserData< TData, TDim > > spData, Grid &grid, const int si)
 
template<typename TFunction , typename TData >
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)
 
void write_item_to_file (VTKFileWriter &File, float data)
 writes data to stream More...
 
void write_item_to_file (VTKFileWriter &File, double data)
 writes data to stream More...
 
void write_item_to_file (VTKFileWriter &File, const ug::MathVector< 1 > &data)
 writes data to stream More...
 
void write_item_to_file (VTKFileWriter &File, const ug::MathVector< 2 > &data)
 writes data to stream More...
 
void write_item_to_file (VTKFileWriter &File, const ug::MathVector< 3 > &data)
 writes data to stream More...
 
void write_item_to_file (VTKFileWriter &File, const ug::MathMatrix< 1, 1 > &data)
 writes data to stream More...
 
void write_item_to_file (VTKFileWriter &File, const ug::MathMatrix< 2, 2 > &data)
 writes data to stream More...
 
void write_item_to_file (VTKFileWriter &File, const ug::MathMatrix< 3, 3 > &data)
 writes data to stream More...
 

Protected Attributes

bool m_bBinary
 print values in binary (base64 encoded way) or plain ascii More...
 
bool m_bSelectAll
 scheduled components to be printed More...
 
bool m_bWriteGrid
 
bool m_bWriteProcRanks
 
bool m_bWriteSubsetIndices
 
std::map< std::string, std::vector< number > > m_mTimestep
 map storing the time points More...
 
std::string m_sComment
 
std::map< std::string, SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > > > m_vMatrixElemData
 
std::map< std::string, SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > > > m_vMatrixNodalData
 scheduled vector data to be printed More...
 
std::map< std::string, SmartPtr< UserData< number, TDim > > > m_vScalarElemData
 
std::map< std::string, SmartPtr< UserData< number, TDim > > > m_vScalarNodalData
 scheduled scalar data to be printed More...
 
std::map< std::string, std::vector< std::string > > m_vSymbFct
 
std::map< std::string, std::vector< std::string > > m_vSymbFctElem
 
std::map< std::string, std::vector< std::string > > m_vSymbFctNodal
 
std::map< std::string, SmartPtr< UserData< MathVector< TDim >, TDim > > > m_vVectorElemData
 
std::map< std::string, SmartPtr< UserData< MathVector< TDim >, TDim > > > m_vVectorNodalData
 scheduled vector data to be printed More...
 

Detailed Description

template<int TDim>
class ug::VTKOutput< TDim >

output writer to the VTK file format

This class can be used to write grid and data associated with the grid to the vtk file format.

The produced files are:

1)) TIME DEPENDENT 1a) In case of functions that are only defined for some subsets

  • filename_pXXXX_sXXXX_tXXXX.vtu (piece of grid for proc, subset, timestep)
  • filename_sXXXX_tXXXX.pvtu (piece of grid for subset, timestep)
  • filename_tXXXX.pvd (piece of grid for timestep)
  • filename.pvd (group of timeseries)

1b) In case that all functions are defined globally

  • filename_p0000_t0000.vtu (piece of grid for proc, timestep)
  • filename_t0000.pvtu (piece of grid for timestep)

2)) TIME INDEPENDENT 2a) In case of functions that are only defined for some subsets

  • filename_pXXXX_sXXXX.vtu (piece of grid for proc, subset)
  • filename_sXXXX.pvtu (piece of grid for subset)
  • filename.pvd (group of subsets)

2b) In case that all functions are defined globally

  • filename_p0000.vtu (piece of grid for proc)
  • filename.pvtu (group of procs)

In case that the simulation is run in serial (or in parallel environment but with only one process), the *.pvtu are not written and the *.vtu file contains all information.

ATTENTION: This class uses heavily the mark-function of the grid. Do not use any member function while having called begin_mark()

Template Parameters
dimworld dimension

Member Typedef Documentation

◆ ComponentsIterator

template<int TDim>
typedef std::map<std::string, std::vector<std::string> >::iterator ug::VTKOutput< TDim >::ComponentsIterator
protected

◆ MatrixDataIterator

template<int TDim>
typedef std::map<std::string, SmartPtr<UserData<MathMatrix<TDim,TDim>, TDim> > >::iterator ug::VTKOutput< TDim >::MatrixDataIterator
protected

◆ ScalarDataIterator

template<int TDim>
typedef std::map<std::string, SmartPtr<UserData<number, TDim> > >::iterator ug::VTKOutput< TDim >::ScalarDataIterator
protected

◆ VectorDataIterator

template<int TDim>
typedef std::map<std::string, SmartPtr<UserData<MathVector<TDim>, TDim> > >::iterator ug::VTKOutput< TDim >::VectorDataIterator
protected

Constructor & Destructor Documentation

◆ VTKOutput()

template<int TDim>
ug::VTKOutput< TDim >::VTKOutput ( )
inline

default constructor

Member Function Documentation

◆ clear_data_selection()

template<int TDim>
void ug::VTKOutput< TDim >::clear_data_selection ( )
inline

◆ clear_selection()

template<int TDim>
void ug::VTKOutput< TDim >::clear_selection ( )
inline

◆ count_piece_sizes()

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::count_piece_sizes ( Grid grid,
const T &  iterContainer,
const SubsetGroup ssGrp,
const int  dim,
int &  numVert,
int &  numElem,
int &  numConn 
)
protected

This function counts the number of vertices, elements and connections for given subsets.

Parameters
[in]udiscrete function
[in]ssGrpsubset group of the piece
[in]dimdimension of subset
[out]numVertnumber of vertices
[out]numElemnumber of elements
[out]numConnnumber of connections

Counts VTK grid elements (points, full-dim. elements and connections) in a given piece of the ug grid.

References ug::Grid::begin_marking(), dim, ug::Grid::end_marking(), ug::SubsetGroup::size(), and UG_THROW.

◆ count_sizes()

template<int TDim>
template<typename TElem , typename T >
void ug::VTKOutput< TDim >::count_sizes ( Grid grid,
const T &  iterContainer,
const int  si,
int &  numVert,
int &  numElem,
int &  numConn 
)
protected

This function counts the number of vertices, elements and connections that are in a subset of the grid (or whole grid, if and only if si < 0).

NOTE: that the number of vertices that are needed to describe the subset may be more than the vertices of the subset. This happens, when elements are part of the subset, but the vertices of these elements are part of another subset. This function counts also those vertices, i.e. it counts the closure of all vertices needed to describe the subset.

Parameters
[in]udiscrete function
[in]sisubset index (si < 0 indicates whole grid)
[out]numVertnumber of vertices needed to compose subset
[out]numElemnumber of elements in the subset
[out]numConnnumber of connections

References ug::GetVertex(), ug::Grid::is_marked(), ug::Grid::mark(), and ug::ROID_OCTAHEDRON.

◆ print() [1/5]

◆ print() [2/5]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print ( const char *  filename,
TFunction &  u 
)
inline

Calls print with makeConsistent enabled.

References ug::VTKOutput< TDim >::print().

◆ print() [3/5]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print ( const char *  filename,
TFunction &  u,
bool  makeConsistent 
)
inline

This function simply calls 'print' using step = -1 to indicate the stationary case. It is intended to write time independent data.

Parameters
[in]filenamefilename for produced files
[in]ugrid function
[in]makeConsistentflag if data shall be made consistent before printing

References ug::VTKOutput< TDim >::print().

◆ print() [4/5]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print ( const char *  filename,
TFunction &  u,
int  step,
number  time 
)
inline

Calls print with makeConsistent enabled.

References ug::VTKOutput< TDim >::print().

◆ print() [5/5]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print ( const char *  filename,
TFunction &  u,
int  step,
number  time,
bool  makeConsistent 
)

This function writes the grid to a *.vtu file. It calls the function print_subset for each subset if necessary and produces a grouping *.pvd file for paraview. If all subsets have to be written, then only a single file is produced containing the whole grid.

Parameters
[in]filenamefilename for produced files
[in]ugrid function
[in]steptime step counter (-1 indicates stationary case)
[in]timetime point corresponding to timestep
[in]makeConsistentflag if data shall be made consistent before printing

VTU output for the data based on the given grid function (u), for all the subsets where this function is defined.

References ug::SubsetGroup::add_all(), PCL_DEBUG_BARRIER_ALL, PROFILE_FUNC, ug::PST_CONSISTENT, UG_CATCH_THROW, and UG_THROW.

Referenced by ug::NestedIterationSolver< TDomain, TAlgebra >::apply(), ug::VTKOutput< TDim >::print(), ug::bridge::Output::SaveDomainToVTK(), ug::GridFunctionDebugWriter< TDomain, TAlgebra >::write_vector_to_vtk(), and ug::bridge::Output::WriteGridFunctionToVTK().

◆ print_subset() [1/2]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print_subset ( const char *  filename,
TFunction &  u,
int  si,
int  step,
number  time,
bool  makeConsistent 
)

This function writes the subset si of the grid (or the whole grid if si < 0) to the file "filename.vtu".

If step >= 0 is passed, this indicates that a time step is written and to the filename a "*_pXXXX*.vtu" is added, indicating the timestep.

If the computation is done in parallel and the number of Processes is greater than one, a "*_pXXXX*.vtu" is added, indicating the process. Then, in addition a "filename*.pvtu" file is written, grouping all *vtu files from different processes.

If only a part of the grid, i.e. a subset, is written, the to the filename a "*_sXXXX*.vtu" is added, to indicate the subset.

Parameters
[in,out]filenamebase name for output file(s)
[in]ugrid function
[in]siSubset (si < 0 indicates whole grid)
[in]stepcounter for timestep (-1 means stationary)
[in]timetime point of timestep
[in]makeConsistentflag if data shall be made consistent before printing

VTU output for the data based on the given grid function (u), for only one given subset si. The file name is composed with the number of the subset.

References ug::Grid::AttachmentAccessor< TElem, TAttachment >::access(), ug::SubsetGroup::add(), ug::Grid::attach_to_vertices(), ug::Grid::detach_from_vertices(), dim, ug::DimensionOfSubset(), ug::IsLittleEndian(), name, ug::Base64FileWriter::normal, pcl::NumProcs(), PCL_DEBUG_BARRIER_ALL, pcl::ProcRank(), PROFILE_FUNC, ug::PST_CONSISTENT, UG_CATCH_THROW, and UG_THROW.

Referenced by ug::VTKOutput< TDim >::print_subset().

◆ print_subset() [2/2]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print_subset ( const char *  filename,
TFunction &  u,
int  si,
int  step = -1,
number  time = 0.0 
)
inline

Calls print_subset with makeConsistent enabled.

References ug::VTKOutput< TDim >::print_subset().

◆ print_subsets() [1/5]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print_subsets ( const char *  filename,
TFunction &  u,
const char *  ssNames 
)
inline

Calls print_subsets with makeConsistent enabled.

References ug::VTKOutput< TDim >::print_subsets().

◆ print_subsets() [2/5]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print_subsets ( const char *  filename,
TFunction &  u,
const char *  ssNames,
bool  makeConsistent 
)
inline

This function simply calls 'print_subsets' using step = -1 to indicate the stationary case. It is intended to write time independent data.

Parameters
[in]filenamefilename for produced files
[in]ugrid function
[in]ssNamessubset names
[in]makeConsistentflag if data shall be made consistent before printing

References ug::VTKOutput< TDim >::print_subsets().

◆ print_subsets() [3/5]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print_subsets ( const char *  filename,
TFunction &  u,
const char *  ssNames,
int  step,
number  time 
)
inline

Calls print_subsets with makeConsistent enabled.

References ug::VTKOutput< TDim >::print_subsets().

◆ print_subsets() [4/5]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print_subsets ( const char *  filename,
TFunction &  u,
const char *  ssNames,
int  step = -1,
number  time = 0.0,
bool  makeConsistent = true 
)

This function writes the given subsets ssGrp of the grid "filename.vtu".

If step >= 0 is passed, this indicates that a time step is written and to the filename a "*_pXXXX*.vtu" is added, indicating the timestep.

If the computation is done in parallel and the number of Processes is greater than one, a "*_pXXXX*.vtu" is added, indicating the process. Then, in addition a "filename*.pvtu" file is written, grouping all *vtu files from different processes.

Parameters
[in,out]filenamebase name for output file(s)
[in]ugrid function
[in]ssNamessubset names
[in]stepcounter for timestep (-1 means stationary)
[in]timetime point of timestep
[in]makeConsistentflag if data shall be made consistent before printing

VTU output for the data based on the given grid function (u), only for the subsets from a given list of subset names.

References ug::SubsetGroup::add(), and ug::TokenizeTrimString().

◆ print_subsets() [5/5]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::print_subsets ( const char *  filename,
TFunction &  u,
const SubsetGroup ssGrp,
int  step = -1,
number  time = 0.0,
bool  makeConsistent = true 
)

This function writes the given subsets ssGrp of the grid "filename.vtu".

If step >= 0 is passed, this indicates that a time step is written and to the filename a "*_pXXXX*.vtu" is added, indicating the timestep.

If the computation is done in parallel and the number of Processes is greater than one, a "*_pXXXX*.vtu" is added, indicating the process. Then, in addition a "filename*.pvtu" file is written, grouping all *vtu files from different processes.

Parameters
[in,out]filenamebase name for output file(s)
[in]ugrid function
[in]ssGrpsubset group
[in]stepcounter for timestep (-1 means stationary)
[in]timetime point of timestep
[in]makeConsistentflag if data shall be made consistent before printing

VTU output for the data based on the given grid function (u), only for the subsets from a given subset group.

References ug::Grid::AttachmentAccessor< TElem, TAttachment >::access(), ug::Grid::attach_to_vertices(), ug::Grid::detach_from_vertices(), dim, ug::DimensionOfSubset(), ug::IsLittleEndian(), name, ug::Base64FileWriter::normal, pcl::NumProcs(), PCL_DEBUG_BARRIER_ALL, pcl::ProcRank(), PROFILE_FUNC, ug::PST_CONSISTENT, ug::SubsetGroup::size(), UG_CATCH_THROW, and UG_THROW.

Referenced by ug::VTKOutput< TDim >::print_subsets().

◆ pvd_filename()

template<int TDim>
void ug::VTKOutput< TDim >::pvd_filename ( std::string &  nameOut,
std::string  nameIn 
)
static

creates the needed pvd file name

References ug::baseName().

◆ pvd_time_filename()

template<int TDim>
void ug::VTKOutput< TDim >::pvd_time_filename ( std::string &  nameOut,
std::string  nameIn,
int  step 
)
static

creates the needed pvd file name to group the time steps

References ug::AppendCounterToString(), and ug::baseName().

◆ pvtu_filename()

template<int TDim>
void ug::VTKOutput< TDim >::pvtu_filename ( std::string &  nameOut,
std::string  nameIn,
int  si,
int  maxSi,
int  step 
)
static

create the needed pvtu file name

References ug::AppendCounterToString(), and ug::baseName().

◆ select() [1/5]

template<int TDim>
void ug::VTKOutput< TDim >::select ( const char *  fctName,
const char *  name 
)

selects a value of a grid function value to be written

This function schedules the component(s) passed by symbolic name(s) of an approximation space to be written to the vtk file under a specified name. The type of the data (nodal/element) will be determined based on the trial space of the components (i.e if continuous available or not)

If more than one component is passed, the data will be interpreted as a vector and dim arguments must be passed.

example: fctName = "p"; name = "pressure" example: fctNames = "u,v,w"; name = "velocity"

Parameters
[in]fctNamesymbolic name of component
[in]namename that will appear in the vtk file for the data

References name, ug::TokenizeString(), and ug::TrimString().

◆ select() [2/5]

template<int TDim>
void ug::VTKOutput< TDim >::select ( const std::vector< std::string > &  vFct,
const char *  name 
)

selects a value of a grid function value to be written

This function schedules the component(s) passed by symbolic name(s) of an approximation space to be written to the vtk file under a specified name. The type of the data (nodal/element) will be determined based on the trial space of the components (i.e if continuous available or not)

If more than one component is passed, the data will be interpreted as a vector and dim arguments must be passed.

example: fctName = "p"; name = "pressure" example: fctNames = "u,v,w"; name = "velocity"

Parameters
[in]fctNamesymbolic name of component
[in]namename that will appear in the vtk file for the data

References name, and UG_THROW.

◆ select() [3/5]

template<int TDim>
void ug::VTKOutput< TDim >::select ( SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > >  spData,
const char *  name 
)

selects a data value to be written

This function schedules a user data to be written to the vtk file under a specified name. The type of the data (nodal/element) will be determined based on the trial space of the components (i.e if continuous available or not)

Parameters
[in]spDatadata to be written
[in]namename that will appear in the vtk file for the data

References name.

◆ select() [4/5]

template<int TDim>
void ug::VTKOutput< TDim >::select ( SmartPtr< UserData< MathVector< TDim >, TDim > >  spData,
const char *  name 
)

selects a data value to be written

This function schedules a user data to be written to the vtk file under a specified name. The type of the data (nodal/element) will be determined based on the trial space of the components (i.e if continuous available or not)

Parameters
[in]spDatadata to be written
[in]namename that will appear in the vtk file for the data

References name.

◆ select() [5/5]

template<int TDim>
void ug::VTKOutput< TDim >::select ( SmartPtr< UserData< number, TDim > >  spData,
const char *  name 
)

selects a data value to be written

This function schedules a user data to be written to the vtk file under a specified name. The type of the data (nodal/element) will be determined based on the trial space of the components (i.e if continuous available or not)

Parameters
[in]spDatadata to be written
[in]namename that will appear in the vtk file for the data

References name.

◆ select_all()

template<int TDim>
void ug::VTKOutput< TDim >::select_all ( bool  bAll)
inline

schedules that all components of the passed discrete functions will be written to file

References ug::VTKOutput< TDim >::m_bSelectAll.

◆ select_element() [1/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_element ( const char *  fctName,
const char *  name 
)

selects a element value of a grid function value to be written

This function schedules the component(s) passed by symbolic name(s) of an approximation space to be written to the vtk file under a specified name. If more than one component is passed, the data will be interpreted as a vector and dim arguments must be passed.

example: fctName = "p"; name = "pressure" example: fctNames = "u,v,w"; name = "velocity"

Parameters
[in]fctNamesymbolic name of component
[in]namename that will appear in the vtk file for the data

References name, ug::TokenizeString(), and ug::TrimString().

◆ select_element() [2/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_element ( const std::vector< std::string > &  vFct,
const char *  name 
)

selects a element value of a grid function value to be written

This function schedules the component(s) passed by symbolic name(s) of an approximation space to be written to the vtk file under a specified name. If more than one component is passed, the data will be interpreted as a vector and dim arguments must be passed.

example: fctName = "p"; name = "pressure" example: fctNames = "u,v,w"; name = "velocity"

Parameters
[in]fctNamesymbolic name of component
[in]namename that will appear in the vtk file for the data

References name, and UG_THROW.

◆ select_element() [3/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_element ( SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > >  spData,
const char *  name 
)

selects a element data value to be written

This function schedules a user data to be written to the vtk file under a specified name.

Parameters
[in]spDatadata to be written
[in]namename that will appear in the vtk file for the data

References name, ug::TrimString(), and UG_THROW.

◆ select_element() [4/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_element ( SmartPtr< UserData< MathVector< TDim >, TDim > >  spData,
const char *  name 
)

selects a element data value to be written

This function schedules a user data to be written to the vtk file under a specified name.

Parameters
[in]spDatadata to be written
[in]namename that will appear in the vtk file for the data

References name, ug::TrimString(), and UG_THROW.

◆ select_element() [5/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_element ( SmartPtr< UserData< number, TDim > >  spData,
const char *  name 
)

selects a element data value to be written

This function schedules a user data to be written to the vtk file under a specified name.

Parameters
[in]spDatadata to be written
[in]namename that will appear in the vtk file for the data

References name, ug::TrimString(), and UG_THROW.

◆ select_nodal() [1/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_nodal ( const char *  fctName,
const char *  name 
)

selects a nodal value of a grid function value to be written

This function schedules the component(s) passed by symbolic name(s) of an approximation space to be written to the vtk file under a specified name. Note, that for the trial space of the component an evaluation of the data at the nodes must be available (continuous). If more than one component is passed, the data will be interpreted as a vector and dim arguments must be passed.

example: fctName = "p"; name = "pressure" example: fctNames = "u,v,w"; name = "velocity"

Parameters
[in]fctNamesymbolic name of component
[in]namename that will appear in the vtk file for the data

References name, ug::TokenizeString(), and ug::TrimString().

◆ select_nodal() [2/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_nodal ( const std::vector< std::string > &  vFct,
const char *  name 
)

selects a nodal value of a grid function value to be written

This function schedules the component(s) passed by symbolic name(s) of an approximation space to be written to the vtk file under a specified name. Note, that for the trial space of the component an evaluation of the data at the nodes must be available (continuous). If more than one component is passed, the data will be interpreted as a vector and dim arguments must be passed.

example: fctName = "p"; name = "pressure" example: fctNames = "u,v,w"; name = "velocity"

Parameters
[in]fctNamesymbolic name of component
[in]namename that will appear in the vtk file for the data

References name, and UG_THROW.

◆ select_nodal() [3/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_nodal ( SmartPtr< UserData< MathMatrix< TDim, TDim >, TDim > >  spData,
const char *  name 
)

selects a nodal data value to be written

This function schedules a user data to be written to the vtk file under a specified name. Note, that for the the data at the nodes must be available (continuous).

Parameters
[in]spDatadata to be written
[in]namename that will appear in the vtk file for the data

References name, ug::TrimString(), and UG_THROW.

◆ select_nodal() [4/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_nodal ( SmartPtr< UserData< MathVector< TDim >, TDim > >  spData,
const char *  name 
)

selects a nodal data value to be written

This function schedules a user data to be written to the vtk file under a specified name. Note, that for the the data at the nodes must be available (continuous).

Parameters
[in]spDatadata to be written
[in]namename that will appear in the vtk file for the data

References name, ug::TrimString(), and UG_THROW.

◆ select_nodal() [5/5]

template<int TDim>
void ug::VTKOutput< TDim >::select_nodal ( SmartPtr< UserData< number, TDim > >  spData,
const char *  name 
)

selects a nodal data value to be written

This function schedules a user data to be written to the vtk file under a specified name. Note, that for the the data at the nodes must be available (continuous).

Parameters
[in]spDatadata to be written
[in]namename that will appear in the vtk file for the data

References name, ug::TrimString(), and UG_THROW.

◆ set_binary()

template<int TDim>
void ug::VTKOutput< TDim >::set_binary ( bool  b)
inline

should values be printed in binary (base64 encoded way ) or plain ascii

References ug::VTKOutput< TDim >::m_bBinary.

◆ set_user_defined_comment()

template<int TDim>
void ug::VTKOutput< TDim >::set_user_defined_comment ( const char *  comment)
inline

◆ set_write_grid()

template<int TDim>
void ug::VTKOutput< TDim >::set_write_grid ( bool  b)
inline

◆ set_write_proc_ranks()

template<int TDim>
void ug::VTKOutput< TDim >::set_write_proc_ranks ( bool  b)
inline

◆ set_write_subset_indices()

template<int TDim>
void ug::VTKOutput< TDim >::set_write_subset_indices ( bool  b)
inline

◆ vtk_name_used()

template<int TDim>
bool ug::VTKOutput< TDim >::vtk_name_used ( const char *  name) const
protected

returns true if name for vtk-component is already used

References name.

◆ vtu_filename()

template<int TDim>
void ug::VTKOutput< TDim >::vtu_filename ( std::string &  nameOut,
std::string  nameIn,
int  rank,
int  si,
int  maxSi,
int  step 
)
static

creates the needed vtu file name

References ug::AppendCounterToString(), ug::baseName(), and pcl::NumProcs().

◆ write_asc_float()

template<int TDim>
VTKFileWriter & ug::VTKOutput< TDim >::write_asc_float ( VTKFileWriter File,
float  data 
)
inlineprotected

prints ascii representation of a float in the Float32 format (a protection against the denormalized floats)

Writes a float in the VTU-compatible ascii format:

◆ write_cell_connectivity() [1/2]

template<int TDim>
template<class TElem , typename T >
void ug::VTKOutput< TDim >::write_cell_connectivity ( VTKFileWriter File,
Grid::VertexAttachmentAccessor< Attachment< int > > &  aaVrtIndex,
Grid grid,
const T &  iterContainer,
const int  si 
)
protected

This method writes the 'connectivity' for each element of a subset. The connectivity are the indices of all vertex the element is formed of.

Parameters
[in]Filefile stream
[in]ugrid function
[in]sisubset index

References ug::Base64FileWriter::base64_binary, ug::Base64FileWriter::normal, ug::ROID_OCTAHEDRON, and ug::ROID_PRISM.

◆ write_cell_connectivity() [2/2]

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::write_cell_connectivity ( VTKFileWriter File,
Grid::VertexAttachmentAccessor< Attachment< int > > &  aaVrtIndex,
Grid grid,
const T &  iterContainer,
const SubsetGroup ssGrp,
const int  dim,
const int  numConn 
)
protected

This method writes the 'connectivity' for each element of a subset. The connectivity are the indices of all vertex the element is formed of.

References ug::Base64FileWriter::base64_binary, dim, ug::Base64FileWriter::normal, ug::SubsetGroup::size(), and UG_THROW.

◆ write_cell_data()

template<int TDim>
template<typename TFunction , typename TData >
void ug::VTKOutput< TDim >::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 
)
protected

◆ write_cell_data_elementwise()

template<int TDim>
template<typename TElem , typename TFunction , typename TData >
void ug::VTKOutput< TDim >::write_cell_data_elementwise ( VTKFileWriter File,
TFunction &  u,
number  time,
SmartPtr< UserData< TData, TDim > >  spData,
Grid grid,
const int  si 
)
protected

◆ write_cell_offsets() [1/2]

template<int TDim>
template<class TElem , typename T >
void ug::VTKOutput< TDim >::write_cell_offsets ( VTKFileWriter File,
const T &  iterContainer,
const int  si,
int &  n 
)
protected

This method writes the 'offset' for each element of a subset.

Parameters
[in]Filefile stream
[in]ugrid function
[in]sisubset index
[in]ncounter for vertices

References ug::Base64FileWriter::base64_binary, ug::Base64FileWriter::normal, and ug::ROID_OCTAHEDRON.

◆ write_cell_offsets() [2/2]

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::write_cell_offsets ( VTKFileWriter File,
const T &  iterContainer,
const SubsetGroup ssGrp,
const int  dim,
const int  numElem 
)
protected

This method writes the 'offset' for each element of a subset.

References ug::Base64FileWriter::base64_binary, dim, ug::Base64FileWriter::normal, ug::SubsetGroup::size(), and UG_THROW.

◆ write_cell_proc_ranks() [1/2]

template<int TDim>
template<class TElem , typename T >
void ug::VTKOutput< TDim >::write_cell_proc_ranks ( VTKFileWriter File,
const T &  iterContainer,
const int  si,
MGSubsetHandler sh 
)
protected

This method writes the proc ranks for each cell.

Parameters
[in]Filefile stream
[in]ugrid function
[in]sisubset index
[in]shmultigrid subset handler

References ug::Base64FileWriter::base64_binary, ug::Base64FileWriter::normal, and pcl::ProcRank().

◆ write_cell_proc_ranks() [2/2]

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::write_cell_proc_ranks ( VTKFileWriter File,
const T &  iterContainer,
const SubsetGroup ssGrp,
const int  dim,
const int  numElem,
MGSubsetHandler sh 
)
protected

This method writes the proc ranks for each cell.

References ug::Base64FileWriter::base64_binary, dim, ug::Base64FileWriter::normal, ug::SubsetGroup::size(), and UG_THROW.

◆ write_cell_subset_names()

template<int TDim>
void ug::VTKOutput< TDim >::write_cell_subset_names ( VTKFileWriter File,
MGSubsetHandler sh 
)
protected

This method writes the subset names.

Parameters
[in]Filefile stream
[in]shmultigrid subset handler

References ug::ISubsetHandler::get_subset_name(), ug::Base64FileWriter::normal, and ug::ISubsetHandler::num_subsets().

◆ write_cell_subsets() [1/2]

template<int TDim>
template<class TElem , typename T >
void ug::VTKOutput< TDim >::write_cell_subsets ( VTKFileWriter File,
const T &  iterContainer,
const int  si,
MGSubsetHandler sh 
)
protected

This method writes the 'region' for each element of a subset. The region is the ug4 subset index for a given element.

Parameters
[in]Filefile stream
[in]ugrid function
[in]sisubset index
[in]shmultigrid subset handler

References ug::Base64FileWriter::base64_binary, ug::ISubsetHandler::get_subset_index(), and ug::Base64FileWriter::normal.

◆ write_cell_subsets() [2/2]

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::write_cell_subsets ( VTKFileWriter File,
const T &  iterContainer,
const SubsetGroup ssGrp,
const int  dim,
const int  numElem,
MGSubsetHandler sh 
)
protected

This method writes the 'region' for each element of a subset. The region is the ug4 subset index for a given element.

References ug::Base64FileWriter::base64_binary, dim, ug::Base64FileWriter::normal, ug::SubsetGroup::size(), and UG_THROW.

◆ write_cell_types() [1/2]

template<int TDim>
template<class TElem , typename T >
void ug::VTKOutput< TDim >::write_cell_types ( VTKFileWriter File,
const T &  iterContainer,
const int  si 
)
protected

This method writes the 'type' for each element of a subset. The type is a vtk-format specified number for a given element type.

Parameters
[in]Filefile stream
[in]ugrid function
[in]sisubset index

References ug::Base64FileWriter::base64_binary, ug::Base64FileWriter::normal, ug::ROID_EDGE, ug::ROID_HEXAHEDRON, ug::ROID_OCTAHEDRON, ug::ROID_PRISM, ug::ROID_PYRAMID, ug::ROID_QUADRILATERAL, ug::ROID_TETRAHEDRON, ug::ROID_TRIANGLE, and UG_THROW.

◆ write_cell_types() [2/2]

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::write_cell_types ( VTKFileWriter File,
const T &  iterContainer,
const SubsetGroup ssGrp,
const int  dim,
const int  numElem 
)
protected

This method writes the 'type' for each element of a subset. The type is a vtk-format specified number for a given element type.

References ug::Base64FileWriter::base64_binary, dim, ug::Base64FileWriter::normal, ug::SubsetGroup::size(), and UG_THROW.

◆ write_cell_values()

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::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 
)
protected

This function writes the values of a function as a <DataArray> field to the vtu file.

Parameters
[in,out]Filefile to write the points
[in]udiscrete function
[in]vFctcomponents to be written
[in]namename to appear for the component
[in]gridGrid
[in]ssGrpsubset
[in]dimdimension of subsets
[in]numVertnumber of vertices

References ug::Base64FileWriter::base64_binary, dim, name, ug::Base64FileWriter::normal, ug::SubsetGroup::size(), and UG_THROW.

◆ write_cell_values_elementwise()

template<int TDim>
template<typename TElem , typename TFunction >
void ug::VTKOutput< TDim >::write_cell_values_elementwise ( VTKFileWriter File,
TFunction &  u,
const std::vector< size_t > &  vFct,
Grid grid,
const int  si 
)
protected

This method writes the cell values of a function to file. If the function group consists of more than on component, than a vector of components is expected.

Parameters
[in]Filefile stream
[in]ugrid function
[in]vFctcomponents to be written
[in]gridGrid
[in]sisubset index

References ug::AveragePositions(), ug::Base64FileWriter::base64_binary, dim, ug::DoFRef(), ug::Provider< TClass >::get(), ug::Base64FileWriter::normal, ug::LocalDoFSet::num_sh(), ug::ROID_OCTAHEDRON, ug::LocalShapeFunctionSet< TDim, TShape, TGrad >::shapes(), UG_CATCH_THROW, and UG_THROW.

◆ write_cell_values_piece()

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::write_cell_values_piece ( VTKFileWriter File,
TFunction &  u,
number  time,
Grid grid,
const SubsetGroup ssGrp,
const int  dim,
const int  numElem 
)
protected

◆ write_cells()

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::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 
)
protected

This function writes the elements that are part of the specified subsets.

Parameters
[in,out]Filefile to write the points
[in]udiscrete function
[in]ssGrpspecified subsets
[in]dimdimension of subsets
[out]numElemnumber of elements
[out]numConnnumber of connections

This function writes the elements that are part of the specified subsets.

References dim, and ug::Base64FileWriter::normal.

◆ write_comment()

template<int TDim>
void ug::VTKOutput< TDim >::write_comment ( VTKFileWriter File)
protected

◆ write_comment_printf()

template<int TDim>
void ug::VTKOutput< TDim >::write_comment_printf ( FILE *  File)
protected

◆ write_empty_grid_piece()

template<int TDim>
void ug::VTKOutput< TDim >::write_empty_grid_piece ( VTKFileWriter File,
bool  binary = true 
)
static

◆ write_grid_piece()

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::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 
)
protected

Composes only the vtu-piece of the description of the unstructured grid.

References dim, ug::Base64FileWriter::normal, and UG_CATCH_THROW.

◆ write_grid_solution_piece()

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::write_grid_solution_piece ( VTKFileWriter File,
Grid::VertexAttachmentAccessor< Attachment< int > > &  aaVrtIndex,
Grid grid,
TFunction &  u,
number  time,
const SubsetGroup ssGrp,
const int  dim 
)
protected

This function writes a piece of the grid to the vtk file. First the geometric data is written (points, cells, connectivity). Then the data associated to the points or cells is written.

Parameters
[in,out]Filefile to write the points
[in]udiscrete function
[in]ssGrpsubsets
[in]dimdimension of subset

This function writes a piece of the grid to the vtk file. First the geometric data is written (points, cells, connectivity). Then the data associated to the points or cells is written.

References ug::LocalFiniteElementProvider::continuous(), dim, ug::Base64FileWriter::normal, and UG_CATCH_THROW.

◆ write_item_to_file() [1/8]

template<int TDim>
void ug::VTKOutput< TDim >::write_item_to_file ( VTKFileWriter File,
const ug::MathMatrix< 1, 1 > &  data 
)
inlineprotected

writes data to stream

The purpose of the function is to convert a double data to binary float (the Float32-format as used in vtk) or ASCII representation. Note that the ASCII output should contain only numbers in the range of the normalized Float32 representation, i.e. with the minimum absolute value 1.1755e-38 = 2^-126. As 'float'-variables can potentially store denormalized numbers (and this is really the case in the practice) with the minimum absolute value 1.4013e-45 = 2^-149 that are not correctly read by some visualization tools (in some operating systems), these function replace those values (lying near to 0) with 0.

◆ write_item_to_file() [2/8]

template<int TDim>
void ug::VTKOutput< TDim >::write_item_to_file ( VTKFileWriter File,
const ug::MathMatrix< 2, 2 > &  data 
)
inlineprotected

writes data to stream

The purpose of the function is to convert a double data to binary float (the Float32-format as used in vtk) or ASCII representation. Note that the ASCII output should contain only numbers in the range of the normalized Float32 representation, i.e. with the minimum absolute value 1.1755e-38 = 2^-126. As 'float'-variables can potentially store denormalized numbers (and this is really the case in the practice) with the minimum absolute value 1.4013e-45 = 2^-149 that are not correctly read by some visualization tools (in some operating systems), these function replace those values (lying near to 0) with 0.

◆ write_item_to_file() [3/8]

template<int TDim>
void ug::VTKOutput< TDim >::write_item_to_file ( VTKFileWriter File,
const ug::MathMatrix< 3, 3 > &  data 
)
inlineprotected

writes data to stream

The purpose of the function is to convert a double data to binary float (the Float32-format as used in vtk) or ASCII representation. Note that the ASCII output should contain only numbers in the range of the normalized Float32 representation, i.e. with the minimum absolute value 1.1755e-38 = 2^-126. As 'float'-variables can potentially store denormalized numbers (and this is really the case in the practice) with the minimum absolute value 1.4013e-45 = 2^-149 that are not correctly read by some visualization tools (in some operating systems), these function replace those values (lying near to 0) with 0.

◆ write_item_to_file() [4/8]

template<int TDim>
void ug::VTKOutput< TDim >::write_item_to_file ( VTKFileWriter File,
const ug::MathVector< 1 > &  data 
)
inlineprotected

writes data to stream

The purpose of the function is to convert a double data to binary float (the Float32-format as used in vtk) or ASCII representation. Note that the ASCII output should contain only numbers in the range of the normalized Float32 representation, i.e. with the minimum absolute value 1.1755e-38 = 2^-126. As 'float'-variables can potentially store denormalized numbers (and this is really the case in the practice) with the minimum absolute value 1.4013e-45 = 2^-149 that are not correctly read by some visualization tools (in some operating systems), these function replace those values (lying near to 0) with 0.

◆ write_item_to_file() [5/8]

template<int TDim>
void ug::VTKOutput< TDim >::write_item_to_file ( VTKFileWriter File,
const ug::MathVector< 2 > &  data 
)
inlineprotected

writes data to stream

The purpose of the function is to convert a double data to binary float (the Float32-format as used in vtk) or ASCII representation. Note that the ASCII output should contain only numbers in the range of the normalized Float32 representation, i.e. with the minimum absolute value 1.1755e-38 = 2^-126. As 'float'-variables can potentially store denormalized numbers (and this is really the case in the practice) with the minimum absolute value 1.4013e-45 = 2^-149 that are not correctly read by some visualization tools (in some operating systems), these function replace those values (lying near to 0) with 0.

◆ write_item_to_file() [6/8]

template<int TDim>
void ug::VTKOutput< TDim >::write_item_to_file ( VTKFileWriter File,
const ug::MathVector< 3 > &  data 
)
inlineprotected

writes data to stream

The purpose of the function is to convert a double data to binary float (the Float32-format as used in vtk) or ASCII representation. Note that the ASCII output should contain only numbers in the range of the normalized Float32 representation, i.e. with the minimum absolute value 1.1755e-38 = 2^-126. As 'float'-variables can potentially store denormalized numbers (and this is really the case in the practice) with the minimum absolute value 1.4013e-45 = 2^-149 that are not correctly read by some visualization tools (in some operating systems), these function replace those values (lying near to 0) with 0.

◆ write_item_to_file() [7/8]

template<int TDim>
void ug::VTKOutput< TDim >::write_item_to_file ( VTKFileWriter File,
double  data 
)
inlineprotected

writes data to stream

The purpose of the function is to convert a double data to binary float (the Float32-format as used in vtk) or ASCII representation. Note that the ASCII output should contain only numbers in the range of the normalized Float32 representation, i.e. with the minimum absolute value 1.1755e-38 = 2^-126. As 'float'-variables can potentially store denormalized numbers (and this is really the case in the practice) with the minimum absolute value 1.4013e-45 = 2^-149 that are not correctly read by some visualization tools (in some operating systems), these function replace those values (lying near to 0) with 0.

References ug::VTKOutput< TDim >::write_item_to_file().

Referenced by ug::VTKOutput< TDim >::write_item_to_file().

◆ write_item_to_file() [8/8]

template<int TDim>
void ug::VTKOutput< TDim >::write_item_to_file ( VTKFileWriter File,
float  data 
)
inlineprotected

writes data to stream

The purpose of the function is to convert a double data to binary float (the Float32-format as used in vtk) or ASCII representation. Note that the ASCII output should contain only numbers in the range of the normalized Float32 representation, i.e. with the minimum absolute value 1.1755e-38 = 2^-126. As 'float'-variables can potentially store denormalized numbers (and this is really the case in the practice) with the minimum absolute value 1.4013e-45 = 2^-149 that are not correctly read by some visualization tools (in some operating systems), these function replace those values (lying near to 0) with 0.

◆ write_nodal_data()

template<int TDim>
template<typename TFunction , typename TData >
void ug::VTKOutput< TDim >::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 
)
protected

◆ write_nodal_data_elementwise()

template<int TDim>
template<typename TElem , typename TFunction , typename TData >
void ug::VTKOutput< TDim >::write_nodal_data_elementwise ( VTKFileWriter File,
TFunction &  u,
number  time,
SmartPtr< UserData< TData, TDim > >  spData,
Grid grid,
const int  si 
)
protected

◆ write_nodal_values() [1/2]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::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 
)
protected

This function writes the values of a function as a <DataArray> field to the vtu file.

Parameters
[in,out]Filefile to write the points
[in]udiscrete function
[in]vFctcomponents to be written
[in]namename to appear for the component
[in]gridGrid
[in]sisubset
[in]dimdimension of subset
[in]numVertnumber of vertices

◆ write_nodal_values() [2/2]

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::write_nodal_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  numVert 
)
protected

This function writes the values of a function as a <DataArray> field to the vtu file.

Parameters
[in,out]Filefile to write the points
[in]udiscrete function
[in]vFctcomponents to be written
[in]namename to appear for the component
[in]gridGrid
[in]ssGrpsubsets
[in]dimdimension of subsets
[in]numVertnumber of vertices

References ug::Base64FileWriter::base64_binary, ug::Grid::begin_marking(), dim, ug::Grid::end_marking(), name, ug::Base64FileWriter::normal, ug::SubsetGroup::size(), and UG_THROW.

◆ write_nodal_values_elementwise()

template<int TDim>
template<typename TElem , typename TFunction >
void ug::VTKOutput< TDim >::write_nodal_values_elementwise ( VTKFileWriter File,
TFunction &  u,
const std::vector< size_t > &  vFct,
Grid grid,
const int  si 
)
protected

This method writes the nodal values of a function to file. If the function group consists of more than on component, than a vector of components is expected. For each function, access to the (unique) nodal value is required.

Parameters
[in]Filefile stream
[in]ugrid function
[in]vFctcomponents to be written
[in]gridGrid
[in]sisubset index

References ug::Base64FileWriter::base64_binary, ug::DoFRef(), ug::GetVertex(), ug::Grid::is_marked(), ug::Grid::mark(), ug::Base64FileWriter::normal, and UG_THROW.

◆ write_nodal_values_piece()

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::write_nodal_values_piece ( VTKFileWriter File,
TFunction &  u,
number  time,
Grid grid,
const SubsetGroup ssGrp,
const int  dim,
const int  numVert 
)
protected

◆ write_points()

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::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 
)
protected

This function writes the vertices of a piece (the specified subsets) of the grid to a vtk file.

Parameters
[in,out]Filefile to write the points
[in]udiscrete function
[in]ssGrpspecified subsets
[in]dimdimension of subsets
[in]numVertnumber of vertices

This function writes the vertices of a piece (the specified subsets) of the grid to a vtk file.

References ug::Base64FileWriter::base64_binary, ug::Grid::begin_marking(), dim, ug::Grid::end_marking(), ug::Base64FileWriter::normal, ug::SubsetGroup::size(), and UG_THROW.

◆ write_points_cells_piece()

template<int TDim>
template<typename T >
void ug::VTKOutput< TDim >::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 
)
protected

This function composes the vtu sections with the coordinates of the grid points and the connectivities (cells), a part of the unstructured grid description in the vtu file.

References dim, and UG_CATCH_THROW.

◆ write_points_elementwise()

template<int TDim>
template<typename TElem , typename T >
void ug::VTKOutput< TDim >::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 
)
protected

This method loops all elements of a subset and writes the vertex coordinates to a file. All vertices of the element are written reguardless if the vertex itself is contained in the subset or another. Written vertices are marked and not written again.

Parameters
[in]Filefile stream
[in]ugrid function
[in]sisubset index
[in]ncounter for vertices

This method loops all elements of a subset and writes the vertex coordinates to a file. All vertices of the element are written reguardless if the vertex itself is contained in the subset or another. Written vertices are marked and not written again.

References ug::Base64FileWriter::base64_binary, ug::GetVertex(), ug::Grid::is_marked(), ug::Grid::mark(), and ug::Base64FileWriter::normal.

◆ write_pvtu()

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::write_pvtu ( TFunction &  u,
const std::string &  filename,
int  si,
int  step,
number  time 
)
protected

◆ write_subset_pvd()

template<int TDim>
void ug::VTKOutput< TDim >::write_subset_pvd ( int  numSubset,
const std::string &  filename,
int  step = -1,
number  time = 0.0 
)
static

writes a grouping *.pvd file, grouping all data from different subsets

References ug::FilenameWithoutPath(), ug::GetLogAssistant(), ug::LogAssistant::is_output_process(), name, pcl::NumProcs(), PCL_DEBUG_BARRIER_ALL, pcl::ProcRank(), and UG_THROW.

◆ write_time_processwise_pvd()

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::write_time_processwise_pvd ( const char *  filename,
TFunction &  u 
)

When a time series has been computed, this function can be used to procduce a grouping *_processwise.pvd file for paraview visualization.

Parameters
[in]filenamefilename used in time series
[in]ugrid function

References ug::FilenameWithoutPath(), ug::GetLogAssistant(), ug::LogAssistant::is_output_process(), name, pcl::NumProcs(), PCL_DEBUG_BARRIER_ALL, and UG_THROW.

◆ write_time_pvd()

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::write_time_pvd ( const char *  filename,
TFunction &  u 
)

When a time series has been computed, this function can be used to procduce a grouping *.pvd file for paraview visualization.

Parameters
[in]filenamefilename used in time series
[in]ugrid function

References ug::FilenameWithoutPath(), ug::GetLogAssistant(), ug::LogAssistant::is_output_process(), name, pcl::NumProcs(), PCL_DEBUG_BARRIER_ALL, and UG_THROW.

◆ write_time_pvd_subset()

template<int TDim>
template<typename TFunction >
void ug::VTKOutput< TDim >::write_time_pvd_subset ( const char *  filename,
TFunction &  u,
int  si 
)

When a time series has been computed, this function can be used to procduce a grouping *.pvd file for paraview visualization.

Parameters
[in]filenamefilename used in time series
[in]ugrid function
[in]sisubset index

References ug::FilenameWithoutPath(), ug::GetLogAssistant(), ug::LogAssistant::is_output_process(), name, pcl::NumProcs(), PCL_DEBUG_BARRIER_ALL, and UG_THROW.

Member Data Documentation

◆ m_bBinary

template<int TDim>
bool ug::VTKOutput< TDim >::m_bBinary
protected

print values in binary (base64 encoded way) or plain ascii

Referenced by ug::VTKOutput< TDim >::set_binary().

◆ m_bSelectAll

template<int TDim>
bool ug::VTKOutput< TDim >::m_bSelectAll
protected

scheduled components to be printed

Referenced by ug::VTKOutput< TDim >::select_all().

◆ m_bWriteGrid

template<int TDim>
bool ug::VTKOutput< TDim >::m_bWriteGrid
protected

◆ m_bWriteProcRanks

template<int TDim>
bool ug::VTKOutput< TDim >::m_bWriteProcRanks
protected

◆ m_bWriteSubsetIndices

template<int TDim>
bool ug::VTKOutput< TDim >::m_bWriteSubsetIndices
protected

◆ m_mTimestep

template<int TDim>
std::map<std::string, std::vector<number> > ug::VTKOutput< TDim >::m_mTimestep
protected

map storing the time points

◆ m_sComment

template<int TDim>
std::string ug::VTKOutput< TDim >::m_sComment
protected

◆ m_vMatrixElemData

template<int TDim>
std::map<std::string, SmartPtr<UserData<MathMatrix<TDim,TDim>, TDim> > > ug::VTKOutput< TDim >::m_vMatrixElemData
protected

◆ m_vMatrixNodalData

template<int TDim>
std::map<std::string, SmartPtr<UserData<MathMatrix<TDim,TDim>, TDim> > > ug::VTKOutput< TDim >::m_vMatrixNodalData
protected

scheduled vector data to be printed

◆ m_vScalarElemData

template<int TDim>
std::map<std::string, SmartPtr<UserData<number, TDim> > > ug::VTKOutput< TDim >::m_vScalarElemData
protected

◆ m_vScalarNodalData

template<int TDim>
std::map<std::string, SmartPtr<UserData<number, TDim> > > ug::VTKOutput< TDim >::m_vScalarNodalData
protected

scheduled scalar data to be printed

Referenced by ug::VTKOutput< TDim >::clear_data_selection().

◆ m_vSymbFct

template<int TDim>
std::map<std::string, std::vector<std::string> > ug::VTKOutput< TDim >::m_vSymbFct
protected

◆ m_vSymbFctElem

template<int TDim>
std::map<std::string, std::vector<std::string> > ug::VTKOutput< TDim >::m_vSymbFctElem
protected

◆ m_vSymbFctNodal

template<int TDim>
std::map<std::string, std::vector<std::string> > ug::VTKOutput< TDim >::m_vSymbFctNodal
protected

◆ m_vVectorElemData

template<int TDim>
std::map<std::string, SmartPtr<UserData<MathVector<TDim>, TDim> > > ug::VTKOutput< TDim >::m_vVectorElemData
protected

◆ m_vVectorNodalData

template<int TDim>
std::map<std::string, SmartPtr<UserData<MathVector<TDim>, TDim> > > ug::VTKOutput< TDim >::m_vVectorNodalData
protected

scheduled vector data to be printed

Referenced by ug::VTKOutput< TDim >::clear_data_selection().


The documentation for this class was generated from the following files: