/*1This file is part of t8code.2t8code is a C library to manage a collection (a forest) of multiple3connected adaptive space-trees of general element classes in parallel.45Copyright (C) 2024 the developers67t8code is free software; you can redistribute it and/or modify8it under the terms of the GNU General Public License as published by9the Free Software Foundation; either version 2 of the License, or10(at your option) any later version.1112t8code is distributed in the hope that it will be useful,13but WITHOUT ANY WARRANTY; without even the implied warranty of14MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the15GNU General Public License for more details.1617You should have received a copy of the GNU General Public License18along with t8code; if not, write to the Free Software Foundation, Inc.,1951 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.20*/2122/**23* \file t8_vtk_writer.h24* This file contains wrappers for the vtk writer functions.25*/2627#ifndef T8_VTK_WRITER_C_INTERFACE_H28#define T8_VTK_WRITER_C_INTERFACE_H2930#include <t8.h>31#include <t8_vtk.h>32#include <t8_forest/t8_forest_types.h>3334#if T8_ENABLE_VTK35#include <vtkUnstructuredGrid.h>36#endif3738T8_EXTERN_C_BEGIN ();3940#if T8_ENABLE_VTK41/**42* Translate a forest into a vtkUnstructuredGrid with respect to the given flags.43* This function uses the vtk library. t8code must be configured with44* "-DT8CODE_ENABLE_VTK=ON" in order to use it.45* \param [in] forest The forest.46* \param[in, out] unstructuredGrid47* \param [in] write_treeid If true, the global tree id is written for each element.48* \param [in] write_mpirank If true, the mpirank is written for each element.49* \param [in] write_level If true, the refinement level is written for each element.50* \param [in] write_element_id If true, the global element id is written for each element.51* \param [in] curved_flag If true, write the elements as curved element types from vtk.52* \param [in] write_ghosts If true, write out ghost elements as well.53* \param [in] num_data Number of user defined double valued data fields to write.54* \param [in] data Array of t8_vtk_data_field_t of length \a num_data55* providing the user defined per element data.56* If scalar and vector fields are used, all scalar fields57* must come first in the array.58*/59void60t8_forest_to_vtkUnstructuredGrid (t8_forest_t forest, vtkSmartPointer<vtkUnstructuredGrid> unstructuredGrid,61const int write_treeid, const int write_mpirank, const int write_level,62const int write_element_id, const int write_ghosts, const int curved_flag,63const int num_data, t8_vtk_data_field_t *data);64#endif6566/** Write the forest in .pvtu file format. Writes one .vtu file per67* process and a meta .pvtu file.68* This function uses the vtk library. t8code must be configured with69* "-DT8CODE_ENABLE_VTK=ON" in order to use it.70* Currently does not support pyramid elements.71* \param [in] forest The forest.72* \param [in] fileprefix The prefix of the output files. The meta file will be named \a fileprefix.pvtu .73* \param [in] write_treeid If true, the global tree id is written for each element.74* \param [in] write_mpirank If true, the mpirank is written for each element.75* \param [in] write_level If true, the refinement level is written for each element.76* \param [in] write_element_id If true, the global element id is written for each element.77* \param [in] curved_flag If true, write the elements as curved element types from vtk.78* \param [in] write_ghosts If true, write out ghost elements as well.79* \param [in] num_data Number of user defined double valued data fields to write.80* \param [in] data Array of t8_vtk_data_field_t of length \a num_data81* providing the user defined per element data.82* If scalar and vector fields are used, all scalar fields83* must come first in the array.84* \return True if successful, false if not (process local).85* \note If t8code was not configured with vtk, use \ref t8_forest_vtk_write_file86*/87int88t8_forest_vtk_write_file_via_API (t8_forest_t forest, const char *fileprefix, const int write_treeid,89const int write_mpirank, const int write_level, const int write_element_id,90const int curved_flag, const int write_ghosts, const int num_data,91t8_vtk_data_field_t *data);9293/** Write the forest in .pvtu file format. Writes one .vtu file per94* process and a meta .pvtu file.95* This function writes ASCII files and can be used when96* t8code is not configured with "-DT8CODE_ENABLE_VTK=ON" and97* \ref t8_forest_vtk_write_file_via_API is not available.98* \param [in] forest The forest.99* \param [in] fileprefix The prefix of the output files.100* \param [in] write_treeid If true, the global tree id is written for each element.101* \param [in] write_mpirank If true, the mpirank is written for each element.102* \param [in] write_level If true, the refinement level is written for each element.103* \param [in] write_element_id If true, the global element id is written for each element.104* \param [in] write_ghosts If true, each process additionally writes its ghost elements.105* For ghost element the treeid is -1.106* \param [in] num_data Number of user defined double valued data fields to write.107* \param [in] data Array of t8_vtk_data_field_t of length \a num_data108* providing the used defined per element data.109* If scalar and vector fields are used, all scalar fields110* must come first in the array.111* \return True if successful, false if not (process local).112*/113int114t8_forest_vtk_write_file (t8_forest_t forest, const char *fileprefix, const int write_treeid, const int write_mpirank,115const int write_level, const int write_element_id, int write_ghosts, const int num_data,116t8_vtk_data_field_t *data);117118/**119* Write the cmesh in .pvtu file format. Writes one .vtu file per120* process and a meta .pvtu file.121* This function uses the vtk library. t8code must be configured with122* "-DT8CODE_ENABLE_VTK=ON" in order to use it.123*124* \param[in] cmesh The cmesh125* \param[in] fileprefix The prefix of the output files126* \param[in] comm The communicator to use127* \return int128* \note If t8code was not configured with vtk, use \ref t8_cmesh_vtk_write_file129*/130int131t8_cmesh_vtk_write_file_via_API (t8_cmesh_t cmesh, const char *fileprefix, sc_MPI_Comm comm);132133/**134* Write the cmesh in .pvtu file format. Writes one .vtu file per135* process and a meta .pvtu file.136* This function writes ASCII files and can be used when137* t8code is not configured with "-DT8CODE_ENABLE_VTK=ON" and138* \ref t8_cmesh_vtk_write_file_via_API is not available.139*140* \param[in] cmesh The cmesh141* \param[in] fileprefix The prefix of the output files142* \return 0 if successful, non-zero otherwise143*/144int145t8_cmesh_vtk_write_file (t8_cmesh_t cmesh, const char *fileprefix);146147T8_EXTERN_C_END ();148#endif /* T8_VTK_WRITER_C_INTERFACE_H */149150151