44 const std::vector<types::global_dof_index> &new_sol_indices)
57 const unsigned int n_independent_variables)
58 :
n_indep(n_independent_variables)
70 std::vector<std::vector<double>>(
n_indep, std::vector<double>(0));
79 const unsigned int n_independent_variables)
81 ,
n_indep(n_independent_variables)
93 std::vector<std::vector<double>>(
n_indep, std::vector<double>(0));
128 dof_handler->get_triangulation().signals.any_change.connect(
161 dof_handler->get_triangulation().signals.any_change.connect(
204 support_point_quadrature,
206 unsigned int n_support_points =
207 dof_handler->get_fe().get_unit_support_points().size();
208 unsigned int n_components =
dof_handler->get_fe(0).n_components();
221 std::vector<unsigned int> current_fe_index(n_components,
224 std::vector<Point<dim>> current_points(n_components,
Point<dim>());
225 for (
unsigned int support_point = 0; support_point < n_support_points;
229 unsigned int component =
230 dof_handler->get_fe().system_to_component_index(support_point).first;
232 current_fe_index[component] = support_point;
245 for (; cell != endc; ++cell)
249 for (
unsigned int support_point = 0; support_point < n_support_points;
253 .system_to_component_index(support_point)
259 location.
distance(current_points[component]))
262 current_points[component] = test_point;
264 current_fe_index[component] = support_point;
270 std::vector<types::global_dof_index> local_dof_indices(
272 std::vector<types::global_dof_index> new_solution_indices;
273 current_cell->get_dof_indices(local_dof_indices);
293 new_solution_indices.reserve(
dof_handler->get_fe(0).n_components());
294 for (
unsigned int component = 0;
298 new_solution_indices.push_back(
299 local_dof_indices[current_fe_index[component]]);
303 new_point_geometry_data(location, current_points, new_solution_indices);
312 data_entry.second.resize(data_entry.second.size() + n_stored);
346 support_point_quadrature,
348 unsigned int n_support_points =
349 dof_handler->get_fe().get_unit_support_points().size();
350 unsigned int n_components =
dof_handler->get_fe(0).n_components();
365 std::vector<typename DoFHandler<dim>::active_cell_iterator> current_cell(
366 locations.size(), cell);
369 std::vector<Point<dim>> temp_points(n_components,
Point<dim>());
370 std::vector<unsigned int> temp_fe_index(n_components, 0);
371 for (
unsigned int support_point = 0; support_point < n_support_points;
375 unsigned int component =
376 dof_handler->get_fe().system_to_component_index(support_point).first;
378 temp_fe_index[component] = support_point;
380 std::vector<std::vector<Point<dim>>> current_points(
381 locations.size(), temp_points);
382 std::vector<std::vector<unsigned int>> current_fe_index(locations.size(),
394 for (; cell != endc; ++cell)
397 for (
unsigned int support_point = 0; support_point < n_support_points;
401 .system_to_component_index(support_point)
406 for (
unsigned int point = 0; point < locations.size(); ++point)
408 if (locations[point].distance(test_point) <
409 locations[point].distance(current_points[point][component]))
412 current_points[point][component] = test_point;
413 current_cell[point] = cell;
414 current_fe_index[point][component] = support_point;
420 std::vector<types::global_dof_index> local_dof_indices(
422 for (
unsigned int point = 0; point < locations.size(); ++point)
424 current_cell[point]->get_dof_indices(local_dof_indices);
425 std::vector<types::global_dof_index> new_solution_indices;
427 new_solution_indices.reserve(
dof_handler->get_fe(0).n_components());
428 for (
unsigned int component = 0;
432 new_solution_indices.push_back(
433 local_dof_indices[current_fe_index[point][component]]);
437 new_point_geometry_data(locations[point],
438 current_points[point],
439 new_solution_indices);
449 data_entry.second.resize(data_entry.second.size() + n_stored);
469 if (mask.represents_the_all_selected_mask() ==
false)
473 std::make_pair(vector_name,
480 std::pair<std::string, std::vector<std::string>> empty_names(
481 vector_name, std::vector<std::string>());
486 std::pair<std::string, std::vector<std::vector<double>>> pair_data;
487 pair_data.first = vector_name;
488 const unsigned int n_stored =
489 (mask.represents_the_all_selected_mask() ==
false ?
490 mask.n_selected_components() :
495 std::vector<std::vector<double>> vector_size(n_datastreams,
496 std::vector<double>(0));
497 pair_data.second = vector_size;
505 const unsigned int n_components)
507 ComponentMask temp_mask(std::vector<bool>(n_components,
true));
515 const std::string &vector_name,
516 const std::vector<std::string> &component_names)
518 typename std::map<std::string, std::vector<std::string>>::iterator names =
523 typename std::map<std::string, ComponentMask>::iterator mask =
526 unsigned int n_stored = mask->second.n_selected_components();
527 Assert(component_names.size() == n_stored,
530 names->second = component_names;
537 const std::vector<std::string> &independent_names)
577template <
typename VectorType>
580 const VectorType &solution)
601 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
602 data_store_field =
data_store.find(vector_name);
606 typename std::map<std::string, ComponentMask>::iterator mask =
610 unsigned int n_stored =
611 mask->second.n_selected_components(
dof_handler->get_fe(0).n_components());
613 typename std::vector<
617 ++point, ++data_store_index)
624 for (
unsigned int store_index = 0, comp = 0;
628 if (mask->second[comp])
630 unsigned int solution_index = point->solution_indices[comp];
632 ->second[data_store_index * n_stored + store_index]
645template <
typename VectorType>
648 const std::vector<std::string> &vector_names,
649 const VectorType &solution,
667 const UpdateFlags update_flags =
672 "The update of normal vectors may not be requested for evaluation of "
673 "data on cells via DataPostprocessor."));
675 unsigned int n_components =
dof_handler->get_fe(0).n_components();
676 unsigned int n_quadrature_points = quadrature.
size();
678 unsigned int n_output_variables = data_postprocessor.
get_names().size();
683 std::vector<typename VectorType::value_type> scalar_solution_values(
684 n_quadrature_points);
685 std::vector<Tensor<1, dim, typename VectorType::value_type>>
686 scalar_solution_gradients(n_quadrature_points);
687 std::vector<Tensor<2, dim, typename VectorType::value_type>>
688 scalar_solution_hessians(n_quadrature_points);
690 std::vector<Vector<typename VectorType::value_type>> vector_solution_values(
693 std::vector<std::vector<Tensor<1, dim, typename VectorType::value_type>>>
694 vector_solution_gradients(
699 std::vector<std::vector<Tensor<2, dim, typename VectorType::value_type>>>
700 vector_solution_hessians(
706 typename std::vector<
711 const auto reference_cell =
712 dof_handler->get_triangulation().get_reference_cells()[0];
714 ++point, ++data_store_index)
717 const Point<dim> requested_location = point->requested_location;
727 std::vector<Vector<double>> computed_quantities(
731 std::vector<Point<dim>> quadrature_points =
733 double distance = cell->diameter();
734 unsigned int selected_point = 0;
735 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
737 if (requested_location.
distance(quadrature_points[q_point]) <
740 selected_point = q_point;
742 requested_location.
distance(quadrature_points[q_point]);
748 if (n_components == 1)
766 std::vector<double>(1, scalar_solution_values[selected_point]);
771 scalar_solution_gradients);
773 std::vector<Tensor<1, dim>>(
774 1, scalar_solution_gradients[selected_point]);
779 scalar_solution_hessians);
781 std::vector<Tensor<2, dim>>(
782 1, scalar_solution_hessians[selected_point]);
787 std::vector<Point<dim>>(1, quadrature_points[selected_point]);
791 computed_quantities);
803 std::copy(vector_solution_values[selected_point].begin(),
804 vector_solution_values[selected_point].end(),
810 vector_solution_gradients);
813 std::copy(vector_solution_gradients[selected_point].begin(),
814 vector_solution_gradients[selected_point].end(),
820 vector_solution_hessians);
823 std::copy(vector_solution_hessians[selected_point].begin(),
824 vector_solution_hessians[selected_point].end(),
829 std::vector<Point<dim>>(1, quadrature_points[selected_point]);
832 computed_quantities);
838 typename std::vector<std::string>::const_iterator name =
839 vector_names.begin();
840 for (; name != vector_names.end(); ++name)
842 typename std::map<std::string,
843 std::vector<std::vector<double>>>::iterator
848 typename std::map<std::string, ComponentMask>::iterator mask =
853 unsigned int n_stored =
854 mask->second.n_selected_components(n_output_variables);
858 for (
unsigned int store_index = 0, comp = 0;
859 comp < n_output_variables;
862 if (mask->second[comp])
865 ->second[data_store_index * n_stored + store_index]
866 .push_back(computed_quantities[0](comp));
877template <
typename VectorType>
880 const std::string &vector_name,
881 const VectorType &solution,
885 std::vector<std::string> vector_names;
886 vector_names.push_back(vector_name);
887 evaluate_field(vector_names, solution, data_postprocessor, quadrature);
893template <
typename VectorType>
896 const std::string &vector_name,
897 const VectorType &solution)
899 using number =
typename VectorType::value_type;
918 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
919 data_store_field =
data_store.find(vector_name);
923 typename std::map<std::string, ComponentMask>::iterator mask =
927 unsigned int n_stored =
928 mask->second.n_selected_components(
dof_handler->get_fe(0).n_components());
930 typename std::vector<
935 ++point, ++data_store_index)
942 point->requested_location,
947 for (
unsigned int store_index = 0, comp = 0; comp < mask->second.size();
950 if (mask->second[comp])
953 ->second[data_store_index * n_stored + store_index]
954 .push_back(
value(comp));
980 const std::vector<double> &indep_values)
993 for (
unsigned int component = 0; component <
n_indep; ++component)
1002 const std::string &base_name,
1003 const std::vector<
Point<dim>> &postprocessor_locations)
1012 std::string filename = base_name +
"_indep.gpl";
1013 std::ofstream to_gnuplot(filename);
1015 to_gnuplot <<
"# Data independent of mesh location\n";
1018 to_gnuplot <<
"# <Key> ";
1024 to_gnuplot <<
"<" << indep_name <<
"> ";
1030 for (
unsigned int component = 0; component <
n_indep; ++component)
1032 to_gnuplot <<
"<Indep_" << component <<
"> ";
1037 for (
unsigned int key = 0; key <
dataset_key.size(); ++key)
1041 for (
unsigned int component = 0; component <
n_indep; ++component)
1058 postprocessor_locations.size() ==
1076 typename std::vector<internal::PointValueHistoryImplementation::
1077 PointGeometryData<dim>>::iterator point =
1079 for (
unsigned int data_store_index = 0;
1081 ++point, ++data_store_index)
1085 std::string filename = base_name +
"_" +
1092 std::ofstream to_gnuplot(filename);
1097 to_gnuplot <<
"# Requested location: " << point->requested_location
1099 to_gnuplot <<
"# DoF_index : Support location (for each component)\n";
1100 for (
unsigned int component = 0;
1101 component <
dof_handler->get_fe(0).n_components();
1104 to_gnuplot <<
"# " << point->solution_indices[component] <<
" : "
1105 << point->support_point_locations[component] <<
'\n';
1109 <<
"# (Original components and locations, may be invalidated by mesh change.)\n";
1111 if (postprocessor_locations.size() != 0)
1113 to_gnuplot <<
"# Postprocessor location: "
1114 << postprocessor_locations[data_store_index];
1116 to_gnuplot <<
" (may be approximate)\n";
1118 to_gnuplot <<
"#\n";
1122 to_gnuplot <<
"# <Key> ";
1128 to_gnuplot <<
"<" << indep_name <<
"> ";
1133 for (
unsigned int component = 0; component <
n_indep; ++component)
1135 to_gnuplot <<
"<Indep_" << component <<
"> ";
1141 typename std::map<std::string, ComponentMask>::iterator mask =
1143 unsigned int n_stored = mask->second.n_selected_components();
1144 std::vector<std::string> names =
1147 if (names.size() > 0)
1151 for (
const auto &name : names)
1153 to_gnuplot <<
"<" << name <<
"> ";
1158 for (
unsigned int component = 0; component < n_stored;
1161 to_gnuplot <<
"<" << data_entry.first <<
"_" << component
1169 for (
unsigned int key = 0; key <
dataset_key.size(); ++key)
1173 for (
unsigned int component = 0; component <
n_indep; ++component)
1180 typename std::map<std::string, ComponentMask>::iterator mask =
1182 unsigned int n_stored = mask->second.n_selected_components();
1184 for (
unsigned int component = 0; component < n_stored;
1189 << (data_entry.second)[data_store_index * n_stored +
1215 typename std::vector<
1220 for (
unsigned int component = 0;
1221 component <
dof_handler->get_fe(0).n_components();
1224 dof_vector(point->solution_indices[component]) = 1;
1234 std::vector<std::vector<
Point<dim>>> &locations)
1240 std::vector<std::vector<Point<dim>>> actual_points;
1241 typename std::vector<
1247 actual_points.push_back(point->support_point_locations);
1249 locations = actual_points;
1263 locations = std::vector<Point<dim>>();
1268 unsigned int n_quadrature_points = quadrature.
size();
1269 std::vector<Point<dim>> evaluation_points;
1272 typename std::vector<
1277 const auto reference_cell =
1278 dof_handler->get_triangulation().get_reference_cells()[0];
1280 ++point, ++data_store_index)
1284 Point<dim> requested_location = point->requested_location;
1291 fe_values.reinit(cell);
1293 evaluation_points = fe_values.get_quadrature_points();
1294 double distance = cell->diameter();
1295 unsigned int selected_point = 0;
1297 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
1299 if (requested_location.
distance(evaluation_points[q_point]) <
1302 selected_point = q_point;
1304 requested_location.
distance(evaluation_points[q_point]);
1308 locations.push_back(evaluation_points[selected_point]);
1317 out <<
"***PointValueHistory status output***\n\n";
1318 out <<
"Closed: " <<
closed <<
'\n';
1319 out <<
"Cleared: " <<
cleared <<
'\n';
1322 out <<
"Geometric Data" <<
'\n';
1324 typename std::vector<
1329 out <<
"No points stored currently\n";
1337 out <<
"# Requested location: " << point->requested_location
1339 out <<
"# DoF_index : Support location (for each component)\n";
1340 for (
unsigned int component = 0;
1341 component <
dof_handler->get_fe(0).n_components();
1344 out << point->solution_indices[component] <<
" : "
1345 << point->support_point_locations[component] <<
'\n';
1352 out <<
"#Cannot access DoF_indices once cleared\n";
1366 out <<
"<" << indep_name <<
"> ";
1373 out <<
"No independent values stored\n";
1379 <<
"Mnemonic: data set size (mask size, n true components) : n data sets\n";
1384 std::string vector_name = data_entry.first;
1385 typename std::map<std::string, ComponentMask>::iterator mask =
1389 typename std::map<std::string, std::vector<std::string>>::iterator
1394 if (data_entry.second.size() != 0)
1396 out << data_entry.first <<
": " << data_entry.second.size() <<
" (";
1397 out << mask->second.size() <<
", "
1398 << mask->second.n_selected_components() <<
") : ";
1399 out << (data_entry.second)[0].size() <<
'\n';
1403 out << data_entry.first <<
": " << data_entry.second.size() <<
" (";
1404 out << mask->second.size() <<
", "
1405 << mask->second.n_selected_components() <<
") : ";
1406 out <<
"No points added" <<
'\n';
1409 if (component_names->second.size() > 0)
1411 for (
const auto &name : component_names->second)
1413 out <<
"<" << name <<
"> ";
1419 out <<
"***end of status output***\n\n";
1444 if ((data_entry.second)[0].size() !=
dataset_key.size())
1471 if (
std::abs(
static_cast<int>((data_entry.second)[0].size()) -
1506#include "numerics/point_value_history.inst"
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
virtual UpdateFlags get_needed_update_flags() const =0
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
virtual void evaluate_scalar_field(const DataPostprocessorInputs::Scalar< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
virtual std::vector< std::string > get_names() const =0
void get_function_values(const ReadVector< Number > &fe_function, std::vector< Number > &values) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
void get_function_hessians(const ReadVector< Number > &fe_function, std::vector< Tensor< 2, spacedim, Number > > &hessians) const
const Point< spacedim > & quadrature_point(const unsigned int q_point) const
void get_function_gradients(const ReadVector< Number > &fe_function, std::vector< Tensor< 1, spacedim, Number > > &gradients) const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
void evaluate_field_at_requested_location(const std::string &name, const VectorType &solution)
void add_field_name(const std::string &vector_name, const ComponentMask &component_mask={})
boost::signals2::connection tria_listener
void add_point(const Point< dim > &location)
std::map< std::string, ComponentMask > component_mask
std::vector< internal::PointValueHistoryImplementation::PointGeometryData< dim > > point_geometry_data
std::map< std::string, std::vector< std::string > > component_names_map
PointValueHistory & operator=(const PointValueHistory &point_value_history)
ObserverPointer< const DoFHandler< dim >, PointValueHistory< dim > > dof_handler
void evaluate_field(const std::string &name, const VectorType &solution)
Vector< double > mark_support_locations()
std::vector< std::string > indep_names
bool triangulation_changed
std::vector< double > dataset_key
void write_gnuplot(const std::string &base_name, const std::vector< Point< dim > > &postprocessor_locations=std::vector< Point< dim > >())
void status(std::ostream &out)
void add_independent_names(const std::vector< std::string > &independent_names)
PointValueHistory(const unsigned int n_independent_variables=0)
void get_support_locations(std::vector< std::vector< Point< dim > > > &locations)
std::map< std::string, std::vector< std::vector< double > > > data_store
void add_component_names(const std::string &vector_name, const std::vector< std::string > &component_names)
void start_new_dataset(const double key)
void add_points(const std::vector< Point< dim > > &locations)
std::vector< std::vector< double > > independent_values
void push_back_independent(const std::vector< double > &independent_values)
void get_postprocessor_locations(const Quadrature< dim > &quadrature, std::vector< Point< dim > > &locations)
void tria_change_listener()
bool deep_check(const bool strict)
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
unsigned int size() const
Point< dim > requested_location
std::vector< types::global_dof_index > solution_indices
std::vector< Point< dim > > support_point_locations
PointGeometryData(const Point< dim > &new_requested_location, const std::vector< Point< dim > > &new_locations, const std::vector< types::global_dof_index > &new_sol_indices)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNoIndependent()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcDataLostSync()
static ::ExceptionBase & ExcDoFHandlerRequired()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDoFHandlerChanged()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
const bool IsBlockVector< VectorType >::value
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
std::vector< Point< spacedim > > evaluation_points
std::vector< double > solution_values
std::vector< Tensor< 2, spacedim > > solution_hessians
std::vector< Tensor< 1, spacedim > > solution_gradients
std::vector< std::vector< Tensor< 2, spacedim > > > solution_hessians
std::vector<::Vector< double > > solution_values
std::vector< std::vector< Tensor< 1, spacedim > > > solution_gradients
static VectorType::value_type get(const VectorType &V, const types::global_dof_index i)