Template Class vtk_writer
Defined in File t8_vtk_writer.hxx
Class Documentation
-
template<typename grid_t>
class vtk_writer A class that controls the writing of vtk files for cmeshes or forests.
- Template Parameters:
grid_t – can be a forest or a cmesh.
Public Functions
-
inline vtk_writer(const bool write_treeid, const bool write_mpirank, const bool write_level, const bool write_element_id, const bool write_ghosts, const bool curved_flag, std::string fileprefix, const int num_data, t8_vtk_data_field_t *data, sc_MPI_Comm comm, const bool merge_points = true)
Construct a new vtk writer object. All parameters are set to false by default. By default no data is used and num_data is set to zero. A default fileprefix is NOT given.
- Parameters:
write_treeid – [in] True, if we want to write the tree id of every element.
write_mpirank – [in] True, if we want to write the mpirankof every element.
write_level – [in] True, if we want to write the level of every element. Uses level 0 if used for a cmesh.
write_element_id – [in] True, if we want to write the element id of every element. Ignored if used for a cmesh.
write_ghosts – [in] True, if we want to write the ghost elements, too.
curved_flag – [in] True, if we want to use quadratic vtk cells. Uses the geometry of the grid to evaluate points between corners.
fileprefix – [in] The prefix of the output-file.
num_data – [in] The number of data-fields to print.
data – [in] The data to use.
comm – [in] The communicator for parallel output.
merge_points – [in] Use the VTK mergePoints feature to merge points shared across elements.
-
inline vtk_writer(std::string fileprefix, sc_MPI_Comm comm)
Construct a new vtk writer object. All parameters are set to false.
- Parameters:
fileprefix – [in]
comm – [in]
-
inline bool write_with_API(const grid_t grid)
A vtk-writer function that uses the vtk API.
- Parameters:
grid – [in] The forest or cmesh that is translated.
- Returns:
true, if writing was successful.
- Returns:
false if writing was not successful.
-
bool write_ASCII(const grid_t grid)
A vtk-writer function that uses the vtk API
- Parameters:
grid – [in] The forest or cmesh that is translated
- Returns:
true
- Returns:
false
-
inline void set_write_treeid(const bool write_treeid)
Set the write treeid flag. Set to true, if you want to write the tree id of every element.
- Parameters:
write_treeid – [in] true or false
-
inline void set_write_mpirank(const bool write_mpirank)
Set the write mpirank flag. Set to true, if you want to write the mpirank of every element.
- Parameters:
write_mpirank – [in] true or false
-
inline void set_write_level(const bool write_level)
Set the write level flag. Set to true, if you want to write the level of every element.
- Parameters:
write_level – [in] true or false
-
inline void set_write_element_id(const bool write_element_id)
Set the write element id flag. Set to true, if you want to write the element id of every element.
- Parameters:
write_element_id – [in] true or false
-
inline void set_write_ghosts(const bool write_ghosts)
Set the write ghosts flag. Set to true, if you want to write the ghost elements, too.
- Parameters:
write_ghosts – [in] true or false
-
inline void set_curved_flag(const bool curved_flag)
Set the curved flag. Set to true, if you want to use quadratic vtk cells. Uses the geometry of the grid to evaluate points between corners.
- Parameters:
curved_flag – [in] true or false
-
inline void set_fileprefix(std::string fileprefix)
Set the fileprefix for the output files.
- Parameters:
fileprefix – [in]
-
bool write_ASCII(const t8_forest_t forest)
Write a forest to a VTK file in ASCII format.
-
bool write_ASCII(const t8_cmesh_t cmesh)
Write a cmesh to a VTK file in ASCII format.