16 #ifndef parsed_tools_data_out_h
17 #define parsed_tools_data_out_h
40 template <
int dim,
int spacedim = dim>
41 class DataOut :
public dealii::ParameterAcceptor
55 const std::string &
base_name =
"solution",
59 const MPI_Comm &
comm = MPI_COMM_WORLD);
69 const std::string &suffix =
"");
80 template <
typename VECTOR>
83 const VECTOR &data_vector,
84 const std::string &desc,
85 const typename dealii::DataOut<dim, spacedim>::DataVectorType &
type =
86 dealii::DataOut<dim, spacedim>::type_automatic);
91 template <
typename VECTOR>
94 const dealii::DataPostprocessor<spacedim> &postproc);
108 dealii::StaticMappingQ1<dim, spacedim>::mapping);
164 std::unique_ptr<dealii::DataOut<dim, spacedim>>
data_out;
166 typename dealii::DataOut<dim, spacedim>::CurvedCellRegion
180 template <
int dim,
int spacedim>
181 template <
typename VECTOR>
184 const VECTOR &data_vector,
185 const std::string &desc,
186 const typename dealii::DataOut<dim, spacedim>::DataVectorType &
type)
188 std::vector<std::string> dd = dealii::Utilities::split_string_list(desc);
189 if (data_out->default_suffix() !=
"")
193 data_out->add_data_vector(data_vector, desc,
type);
197 std::vector<std::string>::iterator sit = dd.begin();
198 std::vector<int> occurrances;
200 for (; sit != dd.end(); ++sit)
201 occurrances.push_back(std::count(dd.begin(), dd.end(), *sit));
203 std::vector<int>::iterator iit = occurrances.begin();
207 dealii::DataComponentInterpretation::DataComponentInterpretation>
208 data_component_interpretation;
210 for (; iit != occurrances.end(); ++iit, ++sit)
213 data_component_interpretation.push_back(
214 dealii::DataComponentInterpretation::
215 component_is_part_of_vector);
217 data_component_interpretation.push_back(
218 dealii::DataComponentInterpretation::component_is_scalar);
221 data_out->add_data_vector(data_vector,
224 data_component_interpretation);
230 template <
int dim,
int spacedim>
231 template <
typename VECTOR>
234 const VECTOR &data_vector,
235 const dealii::DataPostprocessor<spacedim> &postproc)
237 data_out->add_data_vector(data_vector, postproc);
const std::string section_name