16 #include <deal.II/base/quadrature_lib.h> 17 #include <deal.II/base/work_stream.h> 18 #include <deal.II/base/memory_consumption.h> 19 #include <deal.II/lac/vector.h> 20 #include <deal.II/lac/block_vector.h> 21 #include <deal.II/lac/la_vector.h> 22 #include <deal.II/lac/la_parallel_vector.h> 23 #include <deal.II/lac/la_parallel_block_vector.h> 24 #include <deal.II/lac/petsc_vector.h> 25 #include <deal.II/lac/petsc_block_vector.h> 26 #include <deal.II/lac/trilinos_vector.h> 27 #include <deal.II/lac/trilinos_block_vector.h> 28 #include <deal.II/numerics/data_out.h> 29 #include <deal.II/numerics/data_out_dof_data.h> 30 #include <deal.II/dofs/dof_handler.h> 31 #include <deal.II/dofs/dof_accessor.h> 32 #include <deal.II/grid/tria.h> 33 #include <deal.II/grid/tria_iterator.h> 34 #include <deal.II/fe/fe.h> 35 #include <deal.II/fe/fe_dgq.h> 36 #include <deal.II/fe/fe_values.h> 37 #include <deal.II/hp/dof_handler.h> 38 #include <deal.II/hp/fe_values.h> 39 #include <deal.II/fe/mapping_q1.h> 43 DEAL_II_NAMESPACE_OPEN
50 template <
int dim,
int spacedim>
51 ParallelDataBase<dim,spacedim>::
52 ParallelDataBase (
const unsigned int n_datasets,
53 const unsigned int n_subdivisions,
54 const std::vector<unsigned int> &n_postprocessor_outputs,
58 const bool use_face_values)
60 n_datasets (n_datasets),
61 n_subdivisions (n_subdivisions),
62 postprocessed_values (n_postprocessor_outputs.size()),
63 mapping_collection (mapping),
64 finite_elements (finite_elements),
65 update_flags (update_flags)
67 unsigned int n_q_points = 0;
68 if (use_face_values ==
false)
72 n_q_points = quadrature[0].size();
73 x_fe_values.resize(this->finite_elements.size());
74 for (
unsigned int i=0; i<this->finite_elements.size(); ++i)
78 for (
unsigned int j=0; j<i; ++j)
79 if (this->finite_elements[i].
get() ==
80 this->finite_elements[j].get())
82 x_fe_values[i] = x_fe_values[j];
85 if (x_fe_values[i].
get() == 0)
86 x_fe_values[i].reset(new ::hp::FEValues<dim,spacedim>
87 (this->mapping_collection,
88 *this->finite_elements[i],
97 n_q_points = quadrature[0].size();
98 x_fe_face_values.resize(this->finite_elements.size());
99 for (
unsigned int i=0; i<this->finite_elements.size(); ++i)
103 for (
unsigned int j=0; j<i; ++j)
104 if (this->finite_elements[i].
get() ==
105 this->finite_elements[j].get())
107 x_fe_face_values[i] = x_fe_face_values[j];
110 if (x_fe_face_values[i].
get() == 0)
111 x_fe_face_values[i].reset(new ::hp::FEFaceValues<dim,spacedim>
112 (this->mapping_collection,
113 *this->finite_elements[i],
115 this->update_flags));
119 patch_values_scalar.solution_values.resize (n_q_points);
120 patch_values_scalar.solution_gradients.resize (n_q_points);
121 patch_values_scalar.solution_hessians.resize (n_q_points);
122 patch_values_system.solution_values.resize (n_q_points);
123 patch_values_system.solution_gradients.resize (n_q_points);
124 patch_values_system.solution_hessians.resize (n_q_points);
126 for (
unsigned int dataset=0; dataset<n_postprocessor_outputs.size(); ++dataset)
127 if (n_postprocessor_outputs[dataset] != 0)
128 postprocessed_values[dataset]
139 template <
int dim,
int spacedim>
140 ParallelDataBase<dim,spacedim>::
141 ParallelDataBase (
const ParallelDataBase<dim,spacedim> &data)
143 n_datasets (data.n_datasets),
144 n_subdivisions (data.n_subdivisions),
145 patch_values_scalar (data.patch_values_scalar),
146 patch_values_system (data.patch_values_system),
147 postprocessed_values (data.postprocessed_values),
148 mapping_collection (data.mapping_collection),
149 finite_elements (data.finite_elements),
150 update_flags (data.update_flags)
152 if (data.x_fe_values.empty() ==
false)
157 x_fe_values.resize(this->finite_elements.size());
158 for (
unsigned int i=0; i<this->finite_elements.size(); ++i)
162 for (
unsigned int j=0; j<i; ++j)
163 if (this->finite_elements[i].
get() ==
164 this->finite_elements[j].get())
166 x_fe_values[i] = x_fe_values[j];
169 if (x_fe_values[i].
get() == 0)
170 x_fe_values[i].reset(new ::hp::FEValues<dim,spacedim>
171 (this->mapping_collection,
172 *this->finite_elements[i],
174 this->update_flags));
181 x_fe_face_values.resize(this->finite_elements.size());
182 for (
unsigned int i=0; i<this->finite_elements.size(); ++i)
186 for (
unsigned int j=0; j<i; ++j)
187 if (this->finite_elements[i].
get() ==
188 this->finite_elements[j].get())
190 x_fe_face_values[i] = x_fe_face_values[j];
193 if (x_fe_face_values[i].
get() == 0)
194 x_fe_face_values[i].reset(new ::hp::FEFaceValues<dim,spacedim>
195 (this->mapping_collection,
196 *this->finite_elements[i],
198 this->update_flags));
205 template <
int dim,
int spacedim>
206 template <
typename DoFHandlerType>
208 ParallelDataBase<dim,spacedim>::
209 reinit_all_fe_values(std::vector<std_cxx11::shared_ptr<DataEntryBase<DoFHandlerType> > > &dof_data,
210 const typename ::Triangulation<dim,spacedim>::cell_iterator &cell,
211 const unsigned int face)
213 for (
unsigned int dataset=0; dataset<dof_data.size(); ++dataset)
215 bool duplicate =
false;
216 for (
unsigned int j=0; j<dataset; ++j)
217 if (finite_elements[dataset].
get() == finite_elements[j].get())
219 if (duplicate ==
false)
221 typename DoFHandlerType::active_cell_iterator dh_cell(&cell->get_triangulation(),
224 dof_data[dataset]->dof_handler);
225 if (x_fe_values.empty())
229 x_fe_face_values[dataset]->reinit(dh_cell, face);
232 x_fe_values[dataset]->reinit (dh_cell);
235 if (dof_data.empty())
237 if (x_fe_values.empty())
241 x_fe_face_values[0]->reinit(cell, face);
244 x_fe_values[0]->reinit (cell);
250 template <
int dim,
int spacedim>
252 ParallelDataBase<dim,spacedim>::
253 get_present_fe_values(
const unsigned int dataset)
const 256 if (x_fe_values.empty())
257 return x_fe_face_values[dataset]->get_present_fe_values();
259 return x_fe_values[dataset]->get_present_fe_values();
264 template <
int dim,
int spacedim>
266 ParallelDataBase<dim,spacedim>::
271 patch_values_system.solution_gradients.size());
273 patch_values_system.solution_hessians.size());
274 if (patch_values_system.solution_values[0].size() ==
n_components)
276 for (
unsigned int k=0; k<patch_values_system.solution_values.size(); ++k)
278 patch_values_system.solution_values[k].reinit(
n_components);
279 patch_values_system.solution_gradients[k].resize(
n_components);
280 patch_values_system.solution_hessians[k].resize(
n_components);
291 template <
int dim,
int spacedim>
296 patches.push_back (patch);
297 patches.back().patch_index = patches.size()-1;
306 template <
typename DoFHandlerType>
308 (
const DoFHandlerType *dofs,
309 const std::vector<std::string> &names_in,
310 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation)
314 data_component_interpretation (data_component_interpretation),
315 postprocessor(0,
typeid(*this).name()),
316 n_output_variables(names.size())
318 Assert (names.size() == data_component_interpretation.size(),
323 for (
unsigned int i=0; i<names.size(); ++i)
324 Assert (names[i].find_first_not_of(
"abcdefghijklmnopqrstuvwxyz" 325 "ABCDEFGHIJKLMNOPQRSTUVWXYZ" 326 "0123456789_<>()") == std::string::npos,
328 names[i].find_first_not_of(
"abcdefghijklmnopqrstuvwxyz" 329 "ABCDEFGHIJKLMNOPQRSTUVWXYZ" 330 "0123456789_<>()")));
335 template <
typename DoFHandlerType>
337 (
const DoFHandlerType *dofs,
343 postprocessor(data_postprocessor,
typeid(*this).name()),
344 n_output_variables(names.size())
353 for (
unsigned int i=0; i<names.size(); ++i)
354 Assert (names[i].find_first_not_of(
"abcdefghijklmnopqrstuvwxyz" 355 "ABCDEFGHIJKLMNOPQRSTUVWXYZ" 356 "0123456789_<>()") == std::string::npos,
358 names[i].find_first_not_of(
"abcdefghijklmnopqrstuvwxyz" 359 "ABCDEFGHIJKLMNOPQRSTUVWXYZ" 360 "0123456789_<>()")));
365 template <
typename DoFHandlerType>
377 template <
typename DoFHandlerType,
typename VectorType>
387 (
const DoFHandlerType *dofs,
388 const VectorType *data,
389 const std::vector<std::string> &
names,
398 const VectorType *data,
417 std::vector<double> &patch_values)
const;
475 virtual void clear ();
493 template <
typename DoFHandlerType,
typename VectorType>
496 const VectorType *data,
497 const std::vector<std::string> &names,
498 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation)
500 DataEntryBase<DoFHandlerType> (dofs, names, data_component_interpretation),
506 template <
typename DoFHandlerType,
typename VectorType>
509 const VectorType *data,
519 template <
typename VectorType>
521 get_vector_element (
const VectorType &vector,
522 const unsigned int cell_number)
524 return vector[cell_number];
529 get_vector_element (
const IndexSet &is,
530 const unsigned int cell_number)
538 template <
typename DoFHandlerType,
typename VectorType>
543 return get_vector_element(*vector, cell_number);
548 template <
typename DoFHandlerType,
typename VectorType>
565 if (
typeid(
typename VectorType::value_type) ==
typeid(
double))
573 reinterpret_cast<std::vector<::Vector<typename VectorType::value_type>
>&>
574 (patch_values_system));
578 std::vector<::Vector<typename VectorType::value_type> > tmp(patch_values_system.size());
579 for (
unsigned int i = 0; i < patch_values_system.size(); i++)
580 tmp[i].reinit(patch_values_system[i]);
584 for (
unsigned int i = 0; i < patch_values_system.size(); i++)
585 patch_values_system[i] = tmp[i];
591 template <
typename DoFHandlerType,
typename VectorType>
595 std::vector<double> &patch_values)
const 608 if (
typeid(
typename VectorType::value_type) ==
typeid(
double))
616 reinterpret_cast<std::vector<typename VectorType::value_type>&
> 621 std::vector<typename VectorType::value_type> tmp (patch_values.size());
625 for (
unsigned int i = 0; i < tmp.size(); i++)
626 patch_values[i] = tmp[i];
632 template <
typename DoFHandlerType,
typename VectorType>
649 if (
typeid(
typename VectorType::value_type) ==
typeid(
double))
657 reinterpret_cast<std::vector<std::vector<Tensor<1,DoFHandlerType::space_dimension,typename VectorType::value_type>
> >&>
658 (patch_gradients_system));
662 std::vector<std::vector<
Tensor<1,DoFHandlerType::space_dimension,
663 typename VectorType::value_type> > >
664 tmp(patch_gradients_system.size());
665 for (
unsigned int i = 0; i < tmp.size(); i++)
666 tmp[i].resize(patch_gradients_system[i].size());
670 for (
unsigned int i = 0; i < tmp.size(); i++)
671 for (
unsigned int j = 0; j < tmp[i].size(); j++)
672 patch_gradients_system[i][j] = tmp[i][j];
678 template <
typename DoFHandlerType,
typename VectorType>
695 if (
typeid(
typename VectorType::value_type) ==
typeid(
double))
703 reinterpret_cast<std::vector<Tensor<1,DoFHandlerType::space_dimension,typename VectorType::value_type>
>&>
708 std::vector<Tensor<1,DoFHandlerType::space_dimension,typename VectorType::value_type> > tmp;
709 tmp.resize(patch_gradients.size());
713 for (
unsigned int i = 0; i < tmp.size(); i++)
714 patch_gradients[i] = tmp[i];
720 template <
typename DoFHandlerType,
typename VectorType>
737 if (
typeid(
typename VectorType::value_type) ==
typeid(
double))
745 reinterpret_cast<std::vector<std::vector<Tensor<2,DoFHandlerType::space_dimension,typename VectorType::value_type>
> >&>
746 (patch_hessians_system));
750 std::vector<std::vector<
Tensor<2,DoFHandlerType::space_dimension,
751 typename VectorType::value_type> > >
752 tmp(patch_hessians_system.size());
753 for (
unsigned int i = 0; i < tmp.size(); i++)
754 tmp[i].resize(patch_hessians_system[i].size());
758 for (
unsigned int i = 0; i < tmp.size(); i++)
759 for (
unsigned int j = 0; j < tmp[i].size(); j++)
760 patch_hessians_system[i][j] = tmp[i][j];
766 template <
typename DoFHandlerType,
typename VectorType>
783 if (
typeid(
typename VectorType::value_type) ==
typeid(
double))
791 reinterpret_cast<std::vector<
Tensor<2,DoFHandlerType
792 ::space_dimension,typename VectorType::value_type
> >&>
797 std::vector<Tensor<2,DoFHandlerType::space_dimension,typename VectorType::value_type> >
798 tmp(patch_hessians.size());
802 for (
unsigned int i = 0; i < tmp.size(); i++)
803 patch_hessians[i] = tmp[i];
809 template <
typename DoFHandlerType,
typename VectorType>
813 return (
sizeof (vector) +
819 template <
typename DoFHandlerType,
typename VectorType>
824 this->dof_handler = 0;
831 template <
typename DoFHandlerType,
832 int patch_dim,
int patch_space_dim>
835 triangulation(0,typeid(*this).name()),
836 dofs(0,typeid(*this).name())
841 template <
typename DoFHandlerType,
int patch_dim,
int patch_space_dim>
849 template <
typename DoFHandlerType,
int patch_dim,
int patch_space_dim>
854 Assert (dof_data.size() == 0,
856 Assert (cell_data.size() == 0,
860 DoFHandlerType::space_dimension> >
861 (&d.get_triangulation(),
typeid(*this).name());
867 template <
typename DoFHandlerType,
int patch_dim,
int patch_space_dim>
872 Assert (dof_data.size() == 0,
874 Assert (cell_data.size() == 0,
878 DoFHandlerType::space_dimension> >
879 (&tria,
typeid(*this).name());
885 template <
typename DoFHandlerType,
886 int patch_dim,
int patch_space_dim>
887 template <
typename VectorType>
891 const std::string &name,
893 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation)
895 Assert (triangulation != 0,
898 dofs != 0 ? dofs->get_fe().n_components () : 1;
900 std::vector<std::string> names;
903 (vec.size() == triangulation->n_active_cells()))
905 names.resize (1, name);
913 std::ostringstream namebuf;
915 names[i] = name + namebuf.str();
919 add_data_vector (vec, names, type, data_component_interpretation);
924 template <
typename DoFHandlerType,
925 int patch_dim,
int patch_space_dim>
926 template <
typename VectorType>
930 const std::vector<std::string> &names,
932 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation_)
934 Assert (triangulation != 0,
937 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &
938 data_component_interpretation
939 = (data_component_interpretation_.size() != 0
941 data_component_interpretation_
943 std::vector<DataComponentInterpretation::DataComponentInterpretation>
949 if (type == type_automatic)
952 Assert((dofs == 0) || (triangulation->n_active_cells() != dofs->n_dofs()),
953 ExcMessage(
"Unable to determine the type of vector automatically because the number of DoFs " 954 "is equal to the number of cells. Please specify DataVectorType."));
956 if (vec.size() == triangulation->n_active_cells())
957 actual_type = type_cell_data;
959 actual_type = type_dof_data;
965 Assert (vec.size() == triangulation->n_active_cells(),
967 triangulation->n_active_cells()));
968 Assert (names.size() == 1,
975 Assert (vec.size() == dofs->n_dofs(),
978 triangulation->n_active_cells()));
979 Assert (names.size() == dofs->get_fe().n_components(),
981 dofs->get_fe().n_components()));
991 data_component_interpretation);
992 if (actual_type == type_dof_data)
1000 template <
typename DoFHandlerType,
1001 int patch_dim,
int patch_space_dim>
1002 template <
typename VectorType>
1016 Assert (vec.size() == dofs->n_dofs(),
1019 dofs->get_triangulation().n_active_cells()));
1028 template <
typename DoFHandlerType,
1029 int patch_dim,
int patch_space_dim>
1030 template <
typename VectorType>
1034 const VectorType &vec,
1051 template <
typename DoFHandlerType,
1052 int patch_dim,
int patch_space_dim>
1053 template <
typename VectorType>
1057 (
const DoFHandlerType &dof_handler,
1058 const VectorType &data,
1059 const std::string &name,
1060 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation)
1062 const unsigned int n_components = dof_handler.get_fe().n_components ();
1064 std::vector<std::string> names;
1067 names.resize (1, name);
1074 std::ostringstream namebuf;
1075 namebuf <<
'_' << i;
1076 names[i] = name + namebuf.str();
1080 add_data_vector (dof_handler, data, names, data_component_interpretation);
1085 template <
typename DoFHandlerType,
1086 int patch_dim,
int patch_space_dim>
1087 template <
typename VectorType>
1091 (
const DoFHandlerType &dof_handler,
1092 const VectorType &data,
1093 const std::vector<std::string> &names,
1094 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation_)
1099 if (triangulation == 0)
1102 Assert (&dof_handler.get_triangulation() == triangulation,
1103 ExcMessage(
"The triangulation attached to the DoFHandler does not " 1104 "match with the one set previously"));
1106 Assert (data.size() == dof_handler.n_dofs(),
1108 Assert (names.size() == dof_handler.get_fe().n_components(),
1110 dof_handler.get_fe().n_components()));
1112 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &
1113 data_component_interpretation
1114 = (data_component_interpretation_.size() != 0
1116 data_component_interpretation_
1118 std::vector<DataComponentInterpretation::DataComponentInterpretation>
1123 data_component_interpretation);
1129 template <
typename DoFHandlerType,
1130 int patch_dim,
int patch_space_dim>
1133 dof_data.erase (dof_data.begin(), dof_data.end());
1134 cell_data.erase (cell_data.begin(), cell_data.end());
1137 std::vector<Patch> dummy;
1138 patches.swap (dummy);
1143 template <
typename DoFHandlerType,
1144 int patch_dim,
int patch_space_dim>
1149 for (
unsigned int i=0; i<dof_data.size(); ++i)
1150 dof_data[i]->clear ();
1152 for (
unsigned int i=0; i<cell_data.size(); ++i)
1153 cell_data[i]->clear ();
1161 template <
typename DoFHandlerType,
1162 int patch_dim,
int patch_space_dim>
1166 dof_data.erase (dof_data.begin(), dof_data.end());
1167 cell_data.erase (cell_data.begin(), cell_data.end());
1173 std::vector<Patch> dummy;
1174 patches.swap (dummy);
1179 template <
typename DoFHandlerType,
1180 int patch_dim,
int patch_space_dim>
1181 std::vector<std::string>
1185 std::vector<std::string> names;
1189 typename std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > >::const_iterator
1192 for (data_iterator d=dof_data.begin();
1193 d!=dof_data.end(); ++d)
1194 for (
unsigned int i=0; i<(*d)->names.size(); ++i)
1195 names.push_back ((*d)->names[i]);
1196 for (data_iterator d=cell_data.begin(); d!=cell_data.end(); ++d)
1199 names.push_back ((*d)->names[0]);
1207 template <
typename DoFHandlerType,
1208 int patch_dim,
int patch_space_dim>
1209 std::vector<std_cxx11::tuple<unsigned int, unsigned int, std::string> >
1212 std::vector<std_cxx11::tuple<unsigned int, unsigned int, std::string> >
1218 typename std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > >::const_iterator
1221 unsigned int output_component = 0;
1222 for (data_iterator d=dof_data.begin();
1223 d!=dof_data.end(); ++d)
1224 for (
unsigned int i=0; i<(*d)->n_output_variables;
1225 ++i, ++output_component)
1230 if ((*d)->data_component_interpretation[i] ==
1237 Assert (i+patch_space_dim <=
1238 (*d)->n_output_variables,
1241 for (
unsigned int dd=1; dd<patch_space_dim; ++dd)
1242 Assert ((*d)->data_component_interpretation[i+dd]
1254 std::string name = (*d)->names[i];
1255 for (
unsigned int dd=1; dd<patch_space_dim; ++dd)
1256 if (name != (*d)->names[i+dd])
1264 std_cxx11::tuple<unsigned int, unsigned int, std::string>
1265 range (output_component,
1266 output_component+patch_space_dim-1,
1269 ranges.push_back (range);
1275 output_component += patch_space_dim-1;
1276 i += patch_space_dim-1;
1288 unsigned int n_output_components = 0;
1289 for (data_iterator d=dof_data.begin();
1290 d!=dof_data.end(); ++d)
1291 n_output_components += (*d)->n_output_variables;
1292 Assert (output_component == n_output_components,
1301 template <
typename DoFHandlerType,
1302 int patch_dim,
int patch_space_dim>
1303 const std::vector< ::DataOutBase::Patch<patch_dim, patch_space_dim> > &
1311 template <
typename DoFHandlerType,
1312 int patch_dim,
int patch_space_dim>
1314 DoFHandlerType::space_dimension> > >
1317 const unsigned int dhdim = DoFHandlerType::dimension;
1318 const unsigned int dhspacedim = DoFHandlerType::space_dimension;
1319 std::vector<std_cxx11::shared_ptr<::hp::FECollection<dhdim,dhspacedim> > >
1320 finite_elements(this->dof_data.size());
1321 for (
unsigned int i=0; i<this->dof_data.size(); ++i)
1323 Assert (dof_data[i]->dof_handler != 0,
1332 bool duplicate =
false;
1333 for (
unsigned int j=0; j<i; ++j)
1334 if (dof_data[i]->dof_handler == dof_data[j]->dof_handler)
1336 finite_elements[i] = finite_elements[j];
1339 if (duplicate ==
false)
1340 finite_elements[i].reset(new ::hp::FECollection<dhdim,dhspacedim>
1341 (this->dof_data[i]->dof_handler->get_fe()));
1343 if (this->dof_data.empty())
1345 finite_elements.resize(1);
1346 finite_elements[0].reset(new ::hp::FECollection<dhdim,dhspacedim>
1349 return finite_elements;
1354 template <
typename DoFHandlerType,
1355 int patch_dim,
int patch_space_dim>
1367 #include "data_out_dof_data.inst" 1369 DEAL_II_NAMESPACE_CLOSE
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
static ::ExceptionBase & ExcInvalidVectorSize(int arg1, int arg2, int arg3)
void clear_input_data_references()
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInvalidVectorDeclaration(int arg1, std::string arg2)
void attach_triangulation(const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &)
const std::vector< DataComponentInterpretation::DataComponentInterpretation > data_component_interpretation
#define AssertIndexRange(index, range)
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
void attach_dof_handler(const DoFHandlerType &)
static ::ExceptionBase & ExcOldDataStillPresent()
virtual std::vector< std_cxx11::tuple< unsigned int, unsigned int, std::string > > get_vector_data_ranges() const
DataEntryBase(const DoFHandlerType *dofs, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation)
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
virtual std::size_t memory_consumption() const
static ::ExceptionBase & ExcMessage(std::string arg1)
const VectorType * vector
#define Assert(cond, exc)
virtual ~DataOut_DoFData()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
DataEntry(const DoFHandlerType *dofs, const VectorType *data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation)
virtual std::vector< std::string > get_dataset_names() const
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual double get_cell_data_value(const unsigned int cell_number) const
std::size_t memory_consumption() const
virtual void get_function_hessians(const FEValuesBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &fe_patch_values, std::vector< Tensor< 2, DoFHandlerType::space_dimension > > &patch_hessians) const
std::vector< std_cxx11::shared_ptr<::hp::FECollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > get_finite_elements() const
virtual void get_function_gradients(const FEValuesBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &fe_patch_values, std::vector< Tensor< 1, DoFHandlerType::space_dimension > > &patch_gradients) const
virtual std::vector< DataComponentInterpretation::DataComponentInterpretation > get_data_component_interpretation() const
bool is_element(const size_type index) const
static ::ExceptionBase & ExcInvalidCharacter(std::string arg1, size_t arg2)
virtual void get_function_values(const FEValuesBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &fe_patch_values, std::vector< double > &patch_values) const
void clear_data_vectors()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcNoTriangulationSelected()
const std::vector< std::string > names
virtual std::vector< std::string > get_names() const =0
static ::ExceptionBase & ExcNoDoFHandlerSelected()
static ::ExceptionBase & ExcInvalidNumberOfNames(int arg1, int arg2)
static ::ExceptionBase & ExcInternalError()
virtual const std::vector< Patch > & get_patches() const