16 #include <deal.II/numerics/data_out_stack.h> 17 #include <deal.II/base/quadrature_lib.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/dofs/dof_handler.h> 22 #include <deal.II/dofs/dof_accessor.h> 23 #include <deal.II/grid/tria_iterator.h> 24 #include <deal.II/fe/fe.h> 25 #include <deal.II/fe/fe_values.h> 26 #include <deal.II/hp/fe_values.h> 27 #include <deal.II/fe/mapping_q1.h> 31 DEAL_II_NAMESPACE_OPEN
34 template <
int dim,
int spacedim,
typename DoFHandlerType>
44 template <
int dim,
int spacedim,
typename DoFHandlerType>
55 for (
typename std::vector<DataVector>::const_iterator i=
dof_data.begin();
57 Assert (i->data.size() == 0,
59 for (
typename std::vector<DataVector>::const_iterator i=
cell_data.begin();
61 Assert (i->data.size() == 0,
67 template <
int dim,
int spacedim,
typename DoFHandlerType>
78 template <
int dim,
int spacedim,
typename DoFHandlerType>
82 std::vector<std::string> names;
83 names.push_back (name);
88 template <
int dim,
int spacedim,
typename DoFHandlerType>
100 for (std::vector<std::string>::const_iterator name=names.begin(); name!=names.end(); ++name)
102 for (
typename std::vector<DataVector>::const_iterator data_set=
dof_data.begin();
103 data_set!=
dof_data.end(); ++data_set)
104 for (
unsigned int i=0; i<data_set->names.size(); ++i)
107 for (
typename std::vector<DataVector>::const_iterator data_set=
cell_data.begin();
109 for (
unsigned int i=0; i<data_set->names.size(); ++i)
128 template <
int dim,
int spacedim,
typename DoFHandlerType>
129 template <
typename number>
131 const std::string &name)
133 const unsigned int n_components =
dof_handler->get_fe(0).n_components ();
135 std::vector<std::string> names;
139 if ((n_components == 1) ||
142 names.resize (1, name);
148 names.resize (n_components);
149 for (
unsigned int i=0; i<n_components; ++i)
151 std::ostringstream namebuf;
153 names[i] = name + namebuf.str();
161 template <
int dim,
int spacedim,
typename DoFHandlerType>
162 template <
typename number>
164 const std::vector<std::string> &names)
174 (names.size() ==
dof_handler->get_fe(0).n_components())),
177 for (
unsigned int i=0; i<names.size(); ++i)
178 Assert (names[i].find_first_not_of(
"abcdefghijklmnopqrstuvwxyz" 179 "ABCDEFGHIJKLMNOPQRSTUVWXYZ" 180 "0123456789_<>()") == std::string::npos,
182 names[i].find_first_not_of(
"abcdefghijklmnopqrstuvwxyz" 183 "ABCDEFGHIJKLMNOPQRSTUVWXYZ" 184 "0123456789_<>()")));
188 typename std::vector<DataVector>::iterator data_vector=
dof_data.begin();
189 for (; data_vector!=
dof_data.end(); ++data_vector)
190 if (data_vector->names == names)
192 data_vector->data.reinit (vec.
size());
194 data_vector->data.begin());
211 typename std::vector<DataVector>::iterator data_vector=
cell_data.begin();
212 for (; data_vector!=
cell_data.end(); ++data_vector)
213 if (data_vector->names == names)
215 data_vector->data.reinit (vec.
size());
217 data_vector->data.begin());
230 template <
int dim,
int spacedim,
typename DoFHandlerType>
235 unsigned int n_subdivisions = (nnnn_subdivisions != 0)
239 Assert (n_subdivisions >= 1,
246 const unsigned int n_components =
dof_handler->get_fe(0).n_components();
247 const unsigned int n_datasets =
dof_data.size() * n_components +
253 unsigned int n_patches = 0;
254 for (
typename DoFHandlerType::active_cell_iterator
282 const unsigned int n_q_points = patch_points.
size();
283 std::vector<double> patch_values (n_q_points);
284 std::vector<Vector<double> > patch_values_system (n_q_points,
294 default_patch.
data.
reinit (n_datasets, n_q_points*(n_subdivisions+1));
299 typename std::vector< ::DataOutBase::Patch<dim+1,dim+1> >::iterator
301 unsigned int cell_number = 0;
302 for (
typename DoFHandlerType::active_cell_iterator cell=
dof_handler->begin_active();
303 cell !=
dof_handler->end(); ++cell, ++patch, ++cell_number)
305 Assert (cell->is_locally_owned(),
372 x_fe_patch_values.
reinit (cell);
377 for (
unsigned int dataset=0; dataset<
dof_data.size(); ++dataset)
379 if (n_components == 1)
383 for (
unsigned int i=0; i<n_subdivisions+1; ++i)
384 for (
unsigned int q=0; q<n_q_points; ++q)
385 patch->data(dataset,q+n_q_points*i) = patch_values[q];
391 patch_values_system);
392 for (
unsigned int component=0; component<n_components; ++component)
393 for (
unsigned int i=0; i<n_subdivisions+1; ++i)
394 for (
unsigned int q=0; q<n_q_points; ++q)
395 patch->data(dataset*n_components+component,
397 = patch_values_system[q](component);
402 for (
unsigned int dataset=0; dataset<
cell_data.size(); ++dataset)
404 const double value =
cell_data[dataset].data(cell_number);
405 for (
unsigned int q=0; q<n_q_points; ++q)
406 for (
unsigned int i=0; i<n_subdivisions+1; ++i)
407 patch->data(dataset+
dof_data.size()*n_components,
408 q*(n_subdivisions+1)+i) = value;
415 template <
int dim,
int spacedim,
typename DoFHandlerType>
420 for (
typename std::vector<DataVector>::iterator i=
dof_data.begin();
424 for (
typename std::vector<DataVector>::iterator i=
cell_data.begin();
431 template <
int dim,
int spacedim,
typename DoFHandlerType>
446 template <
int dim,
int spacedim,
typename DoFHandlerType>
447 const std::vector< ::DataOutBase::Patch<dim+1,dim+1> > &
455 template <
int dim,
int spacedim,
typename DoFHandlerType>
458 std::vector<std::string> names;
459 for (
typename std::vector<DataVector>::const_iterator dataset=
dof_data.begin();
461 names.insert (names.end(), dataset->names.begin(), dataset->names.end());
462 for (
typename std::vector<DataVector>::const_iterator dataset=
cell_data.begin();
464 names.insert (names.end(), dataset->names.begin(), dataset->names.end());
472 #include "data_out_stack.inst" 475 DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNameAlreadyUsed(std::string arg1)
std::vector< ::DataOutBase::Patch< dim+1, dim+1 > > patches
std::size_t memory_consumption() const
SmartPointer< const DoFHandlerType, DataOutStack< dim, spacedim, DoFHandlerType > > dof_handler
std::vector< DataVector > cell_data
static ::ExceptionBase & ExcDataAlreadyAdded()
void attach_dof_handler(const DoFHandlerType &dof_handler)
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
std::vector< DataVector > dof_data
unsigned int default_subdivisions
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
static ::ExceptionBase & ExcDataNotCleared()
const ::FEValues< dim, spacedim > & get_present_fe_values() const
std::size_t memory_consumption() const
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, lda > > cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
std::vector< std::string > names
static ::ExceptionBase & ExcNoDoFHandlerSelected()
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
void declare_data_vector(const std::string &name, const VectorType vector_type)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcInvalidNumberOfNames(int arg1, int arg2)
static ::ExceptionBase & ExcInvalidCharacter(std::string arg1, size_t arg2)
unsigned int n_subdivisions
static ::ExceptionBase & ExcVectorNotDeclared(std::string arg1)
unsigned int size() const
virtual std::vector< std::string > get_dataset_names() const
void build_patches(const unsigned int n_subdivisions=0)
virtual const std::vector< ::DataOutBase::Patch< dim+1, dim+1 > > & get_patches() const
void validate_dataset_names() const
static ::ExceptionBase & ExcNotImplemented()
void finish_parameter_value()
void new_parameter_value(const double parameter_value, const double parameter_step)
void add_data_vector(const Vector< number > &vec, const std::string &name)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcInternalError()