41 namespace PointValueHistoryImplementation
47 const std::vector<types::global_dof_index> &new_sol_indices)
49 requested_location = new_requested_location;
50 support_point_locations = new_locations;
51 solution_indices = new_sol_indices;
60 const unsigned int n_independent_variables)
61 : n_indep(n_independent_variables)
73 std::vector<std::vector<double>>(
n_indep, std::vector<double>(0));
82 const unsigned int n_independent_variables)
83 : dof_handler(&dof_handler)
84 , n_indep(n_independent_variables)
96 std::vector<std::vector<double>>(
n_indep, std::vector<double>(0));
117 closed = point_value_history.
closed;
118 cleared = point_value_history.
cleared;
124 n_indep = point_value_history.
n_indep;
128 if (have_dof_handler)
131 dof_handler->get_triangulation().signals.any_change.connect(
132 [
this]() { this->tria_change_listener(); });
150 closed = point_value_history.
closed;
151 cleared = point_value_history.
cleared;
157 n_indep = point_value_history.
n_indep;
161 if (have_dof_handler)
164 dof_handler->get_triangulation().signals.any_change.connect(
165 [
this]() { this->tria_change_listener(); });
176 if (have_dof_handler)
178 tria_listener.disconnect();
192 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
193 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
205 dof_handler->get_fe().get_unit_support_points());
207 support_point_quadrature,
209 unsigned int n_support_points =
210 dof_handler->get_fe().get_unit_support_points().size();
211 unsigned int n_components = dof_handler->get_fe(0).n_components();
216 dof_handler->begin_active();
224 std::vector<unsigned int> current_fe_index(n_components,
227 std::vector<Point<dim>> current_points(n_components,
Point<dim>());
228 for (
unsigned int support_point = 0; support_point < n_support_points;
233 unsigned int component =
234 dof_handler->get_fe().system_to_component_index(support_point).first;
236 current_fe_index[component] = support_point;
249 for (; cell != endc; ++cell)
253 for (
unsigned int support_point = 0; support_point < n_support_points;
256 unsigned int component = dof_handler->get_fe()
257 .system_to_component_index(support_point)
263 location.
distance(current_points[component]))
266 current_points[component] = test_point;
268 current_fe_index[component] = support_point;
274 std::vector<types::global_dof_index> local_dof_indices(
275 dof_handler->get_fe().n_dofs_per_cell());
276 std::vector<types::global_dof_index> new_solution_indices;
277 current_cell->get_dof_indices(local_dof_indices);
297 for (
unsigned int component = 0;
298 component < dof_handler->get_fe(0).n_components();
301 new_solution_indices.push_back(
302 local_dof_indices[current_fe_index[component]]);
306 new_point_geometry_data(location, current_points, new_solution_indices);
307 point_geometry_data.push_back(new_point_geometry_data);
309 for (
auto &data_entry : data_store)
313 (component_mask.find(data_entry.first))->
second;
315 data_entry.second.resize(data_entry.second.size() + n_stored);
333 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
334 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
347 dof_handler->get_fe().get_unit_support_points());
349 support_point_quadrature,
351 unsigned int n_support_points =
352 dof_handler->get_fe().get_unit_support_points().size();
353 unsigned int n_components = dof_handler->get_fe(0).n_components();
358 dof_handler->begin_active();
368 std::vector<typename DoFHandler<dim>::active_cell_iterator> current_cell(
369 locations.size(), cell);
372 std::vector<Point<dim>> temp_points(n_components,
Point<dim>());
373 std::vector<unsigned int> temp_fe_index(n_components, 0);
374 for (
unsigned int support_point = 0; support_point < n_support_points;
379 unsigned int component =
380 dof_handler->get_fe().system_to_component_index(support_point).first;
382 temp_fe_index[component] = support_point;
384 std::vector<std::vector<Point<dim>>> current_points(
385 locations.size(), temp_points);
386 std::vector<std::vector<unsigned int>> current_fe_index(locations.size(),
398 for (; cell != endc; ++cell)
401 for (
unsigned int support_point = 0; support_point < n_support_points;
404 unsigned int component = dof_handler->get_fe()
405 .system_to_component_index(support_point)
410 for (
unsigned int point = 0; point < locations.size(); ++point)
412 if (locations[point].distance(test_point) <
413 locations[point].distance(current_points[point][component]))
416 current_points[point][component] = test_point;
417 current_cell[point] = cell;
418 current_fe_index[point][component] = support_point;
424 std::vector<types::global_dof_index> local_dof_indices(
425 dof_handler->get_fe().n_dofs_per_cell());
426 for (
unsigned int point = 0; point < locations.size(); ++point)
428 current_cell[point]->get_dof_indices(local_dof_indices);
429 std::vector<types::global_dof_index> new_solution_indices;
431 for (
unsigned int component = 0;
432 component < dof_handler->get_fe(0).n_components();
435 new_solution_indices.push_back(
436 local_dof_indices[current_fe_index[point][component]]);
440 new_point_geometry_data(locations[point],
441 current_points[point],
442 new_solution_indices);
444 point_geometry_data.push_back(new_point_geometry_data);
446 for (
auto &data_entry : data_store)
450 (component_mask.find(data_entry.first))->
second;
452 data_entry.second.resize(data_entry.second.size() + n_stored);
468 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
469 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
472 if (mask.represents_the_all_selected_mask() ==
false)
473 component_mask.insert(std::make_pair(vector_name, mask));
475 component_mask.insert(
476 std::make_pair(vector_name,
478 dof_handler->get_fe(0).n_components(),
true))));
483 std::pair<std::string, std::vector<std::string>> empty_names(
484 vector_name, std::vector<std::string>());
485 component_names_map.insert(empty_names);
489 std::pair<std::string, std::vector<std::vector<double>>> pair_data;
490 pair_data.first = vector_name;
491 const unsigned int n_stored =
492 (mask.represents_the_all_selected_mask() ==
false ?
493 mask.n_selected_components() :
494 dof_handler->get_fe(0).n_components());
497 point_geometry_data.size() * n_stored;
498 std::vector<std::vector<double>> vector_size(n_datastreams,
499 std::vector<double>(0));
500 pair_data.second = vector_size;
501 data_store.insert(pair_data);
508 const unsigned int n_components)
510 std::vector<bool> temp_mask(n_components,
true);
511 add_field_name(vector_name, temp_mask);
518 const std::string & vector_name,
519 const std::vector<std::string> &component_names)
521 typename std::map<std::string, std::vector<std::string>>::iterator names =
522 component_names_map.find(vector_name);
523 Assert(names != component_names_map.end(),
526 typename std::map<std::string, ComponentMask>::iterator mask =
527 component_mask.find(vector_name);
528 Assert(mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
529 unsigned int n_stored = mask->second.n_selected_components();
531 Assert(component_names.size() == n_stored,
534 names->second = component_names;
541 const std::vector<std::string> &independent_names)
543 Assert(independent_names.size() == n_indep,
546 indep_names = independent_names;
564 dof_handler =
nullptr;
565 have_dof_handler =
false;
581template <
typename VectorType>
584 const VectorType & solution)
590 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
591 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
597 static_cast<int>(independent_values[0].size())) < 2,
605 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
606 data_store_field = data_store.find(vector_name);
607 Assert(data_store_field != data_store.end(),
610 typename std::map<std::string, ComponentMask>::iterator mask =
611 component_mask.find(vector_name);
612 Assert(mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
614 unsigned int n_stored =
615 mask->second.n_selected_components(dof_handler->get_fe(0).n_components());
617 typename std::vector<
619 point = point_geometry_data.begin();
620 for (
unsigned int data_store_index = 0; point != point_geometry_data.end();
621 ++point, ++data_store_index)
628 for (
unsigned int store_index = 0, comp = 0;
629 comp < dof_handler->get_fe(0).n_components();
632 if (mask->second[comp])
634 unsigned int solution_index = point->solution_indices[comp];
636 ->second[data_store_index * n_stored + store_index]
649template <
typename VectorType>
652 const std::vector<std::string> &vector_names,
653 const VectorType & solution,
661 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
666 static_cast<int>(independent_values[0].size())) < 2,
676 "The update of normal vectors may not be requested for evaluation of "
677 "data on cells via DataPostprocessor."));
678 FEValues<dim> fe_values(dof_handler->get_fe(), quadrature, update_flags);
679 unsigned int n_components = dof_handler->get_fe(0).n_components();
680 unsigned int n_quadrature_points = quadrature.
size();
682 unsigned int n_output_variables = data_postprocessor.
get_names().size();
687 std::vector<typename VectorType::value_type> scalar_solution_values(
688 n_quadrature_points);
689 std::vector<Tensor<1, dim, typename VectorType::value_type>>
690 scalar_solution_gradients(n_quadrature_points);
691 std::vector<Tensor<2, dim, typename VectorType::value_type>>
692 scalar_solution_hessians(n_quadrature_points);
694 std::vector<Vector<typename VectorType::value_type>> vector_solution_values(
697 std::vector<std::vector<Tensor<1, dim, typename VectorType::value_type>>>
698 vector_solution_gradients(
703 std::vector<std::vector<Tensor<2, dim, typename VectorType::value_type>>>
704 vector_solution_hessians(
710 typename std::vector<
712 point = point_geometry_data.begin();
713 for (
unsigned int data_store_index = 0; point != point_geometry_data.end();
714 ++point, ++data_store_index)
717 const Point<dim> requested_location = point->requested_location;
726 std::vector<Vector<double>> computed_quantities(
730 std::vector<Point<dim>> quadrature_points =
732 double distance = cell->diameter();
733 unsigned int selected_point = 0;
734 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
736 if (requested_location.
distance(quadrature_points[q_point]) <
739 selected_point = q_point;
741 requested_location.
distance(quadrature_points[q_point]);
747 if (n_components == 1)
765 std::vector<double>(1, scalar_solution_values[selected_point]);
770 scalar_solution_gradients);
772 std::vector<Tensor<1, dim>>(
773 1, scalar_solution_gradients[selected_point]);
778 scalar_solution_hessians);
780 std::vector<Tensor<2, dim>>(
781 1, scalar_solution_hessians[selected_point]);
786 std::vector<Point<dim>>(1, quadrature_points[selected_point]);
790 computed_quantities);
802 std::copy(vector_solution_values[selected_point].begin(),
803 vector_solution_values[selected_point].end(),
809 vector_solution_gradients);
812 std::copy(vector_solution_gradients[selected_point].begin(),
813 vector_solution_gradients[selected_point].end(),
819 vector_solution_hessians);
822 std::copy(vector_solution_hessians[selected_point].begin(),
823 vector_solution_hessians[selected_point].end(),
828 std::vector<Point<dim>>(1, quadrature_points[selected_point]);
831 computed_quantities);
837 typename std::vector<std::string>::const_iterator name =
838 vector_names.begin();
839 for (; name != vector_names.end(); ++name)
841 typename std::map<std::string,
842 std::vector<std::vector<double>>>::iterator
843 data_store_field = data_store.find(*name);
844 Assert(data_store_field != data_store.end(),
847 typename std::map<std::string, ComponentMask>::iterator mask =
848 component_mask.find(*name);
849 Assert(mask != component_mask.end(),
852 unsigned int n_stored =
853 mask->second.n_selected_components(n_output_variables);
857 for (
unsigned int store_index = 0, comp = 0;
858 comp < n_output_variables;
861 if (mask->second[comp])
864 ->second[data_store_index * n_stored + store_index]
865 .push_back(computed_quantities[0](comp));
876template <
typename VectorType>
879 const std::string & vector_name,
880 const VectorType & solution,
884 std::vector<std::string> vector_names;
885 vector_names.push_back(vector_name);
886 evaluate_field(vector_names, solution, data_postprocessor, quadrature);
892template <
typename VectorType>
895 const std::string &vector_name,
896 const VectorType & solution)
898 using number =
typename VectorType::value_type;
903 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
909 static_cast<int>(independent_values[0].size())) < 2,
917 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
918 data_store_field = data_store.find(vector_name);
919 Assert(data_store_field != data_store.end(),
922 typename std::map<std::string, ComponentMask>::iterator mask =
923 component_mask.find(vector_name);
924 Assert(mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
926 unsigned int n_stored =
927 mask->second.n_selected_components(dof_handler->get_fe(0).n_components());
929 typename std::vector<
931 point = point_geometry_data.begin();
933 for (
unsigned int data_store_index = 0; point != point_geometry_data.end();
934 ++point, ++data_store_index)
941 point->requested_location,
946 for (
unsigned int store_index = 0, comp = 0; comp < mask->second.size();
949 if (mask->second[comp])
952 ->second[data_store_index * n_stored + store_index]
953 .push_back(value(comp));
969 Assert(deep_check(
false), ExcDataLostSync());
971 dataset_key.push_back(key);
979 const std::vector<double> &indep_values)
985 Assert(indep_values.size() == n_indep,
987 Assert(n_indep != 0, ExcNoIndependent());
989 static_cast<int>(independent_values[0].size())) < 2,
992 for (
unsigned int component = 0; component < n_indep; ++component)
993 independent_values[component].push_back(indep_values[component]);
1001 const std::string & base_name,
1002 const std::vector<
Point<dim>> &postprocessor_locations)
1011 std::string filename = base_name +
"_indep.gpl";
1012 std::ofstream to_gnuplot(filename.c_str());
1014 to_gnuplot <<
"# Data independent of mesh location\n";
1017 to_gnuplot <<
"# <Key> ";
1019 if (indep_names.size() > 0)
1021 for (
const auto &indep_name : indep_names)
1023 to_gnuplot <<
"<" << indep_name <<
"> ";
1029 for (
unsigned int component = 0; component < n_indep; ++component)
1031 to_gnuplot <<
"<Indep_" << component <<
"> ";
1036 for (
unsigned int key = 0; key < dataset_key.size(); ++key)
1038 to_gnuplot << dataset_key[key];
1040 for (
unsigned int component = 0; component < n_indep; ++component)
1042 to_gnuplot <<
" " << independent_values[component][key];
1053 if (have_dof_handler)
1055 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1056 AssertThrow(postprocessor_locations.size() == 0 ||
1057 postprocessor_locations.size() ==
1058 point_geometry_data.size(),
1060 point_geometry_data.size()));
1075 typename std::vector<internal::PointValueHistoryImplementation::
1076 PointGeometryData<dim>>::iterator point =
1077 point_geometry_data.begin();
1078 for (
unsigned int data_store_index = 0;
1079 point != point_geometry_data.end();
1080 ++point, ++data_store_index)
1084 std::string filename = base_name +
"_" +
1091 std::ofstream to_gnuplot(filename.c_str());
1096 to_gnuplot <<
"# Requested location: " << point->requested_location
1098 to_gnuplot <<
"# DoF_index : Support location (for each component)\n";
1099 for (
unsigned int component = 0;
1100 component < dof_handler->get_fe(0).n_components();
1103 to_gnuplot <<
"# " << point->solution_indices[component] <<
" : "
1104 << point->support_point_locations[component] <<
'\n';
1106 if (triangulation_changed)
1108 <<
"# (Original components and locations, may be invalidated by mesh change.)\n";
1110 if (postprocessor_locations.size() != 0)
1112 to_gnuplot <<
"# Postprocessor location: "
1113 << postprocessor_locations[data_store_index];
1114 if (triangulation_changed)
1115 to_gnuplot <<
" (may be approximate)\n";
1117 to_gnuplot <<
"#\n";
1121 to_gnuplot <<
"# <Key> ";
1123 if (indep_names.size() > 0)
1125 for (
const auto &indep_name : indep_names)
1127 to_gnuplot <<
"<" << indep_name <<
"> ";
1132 for (
unsigned int component = 0; component < n_indep; ++component)
1134 to_gnuplot <<
"<Indep_" << component <<
"> ";
1138 for (
const auto &data_entry : data_store)
1140 typename std::map<std::string, ComponentMask>::iterator mask =
1141 component_mask.find(data_entry.first);
1142 unsigned int n_stored = mask->second.n_selected_components();
1143 std::vector<std::string> names =
1144 (component_names_map.find(data_entry.first))->
second;
1146 if (names.size() > 0)
1150 for (
const auto &name : names)
1152 to_gnuplot <<
"<" << name <<
"> ";
1157 for (
unsigned int component = 0; component < n_stored;
1160 to_gnuplot <<
"<" << data_entry.first <<
"_" << component
1168 for (
unsigned int key = 0; key < dataset_key.size(); ++key)
1170 to_gnuplot << dataset_key[key];
1172 for (
unsigned int component = 0; component < n_indep; ++component)
1174 to_gnuplot <<
" " << independent_values[component][key];
1177 for (
const auto &data_entry : data_store)
1179 typename std::map<std::string, ComponentMask>::iterator mask =
1180 component_mask.find(data_entry.first);
1181 unsigned int n_stored = mask->second.n_selected_components();
1183 for (
unsigned int component = 0; component < n_stored;
1188 << (data_entry.second)[data_store_index * n_stored +
1209 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1210 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
1214 typename std::vector<
1216 point = point_geometry_data.begin();
1217 for (; point != point_geometry_data.end(); ++point)
1219 for (
unsigned int component = 0;
1220 component < dof_handler->get_fe(0).n_components();
1223 dof_vector(point->solution_indices[component]) = 1;
1233 std::vector<std::vector<
Point<dim>>> &locations)
1236 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1237 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
1239 std::vector<std::vector<Point<dim>>> actual_points;
1240 typename std::vector<
1242 point = point_geometry_data.begin();
1244 for (; point != point_geometry_data.end(); ++point)
1246 actual_points.push_back(point->support_point_locations);
1248 locations = actual_points;
1260 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1262 locations = std::vector<Point<dim>>();
1267 unsigned int n_quadrature_points = quadrature.
size();
1268 std::vector<Point<dim>> evaluation_points;
1271 typename std::vector<
1273 point = point_geometry_data.begin();
1274 for (
unsigned int data_store_index = 0; point != point_geometry_data.end();
1275 ++point, ++data_store_index)
1279 Point<dim> requested_location = point->requested_location;
1285 fe_values.reinit(cell);
1287 evaluation_points = fe_values.get_quadrature_points();
1288 double distance = cell->diameter();
1289 unsigned int selected_point = 0;
1291 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
1293 if (requested_location.
distance(evaluation_points[q_point]) <
1296 selected_point = q_point;
1298 requested_location.
distance(evaluation_points[q_point]);
1302 locations.push_back(evaluation_points[selected_point]);
1311 out <<
"***PointValueHistory status output***\n\n";
1312 out <<
"Closed: " << closed <<
'\n';
1313 out <<
"Cleared: " << cleared <<
'\n';
1314 out <<
"Triangulation_changed: " << triangulation_changed <<
'\n';
1315 out <<
"Have_dof_handler: " << have_dof_handler <<
'\n';
1316 out <<
"Geometric Data" <<
'\n';
1318 typename std::vector<
1320 point = point_geometry_data.begin();
1321 if (point == point_geometry_data.end())
1323 out <<
"No points stored currently\n";
1329 for (; point != point_geometry_data.end(); ++point)
1331 out <<
"# Requested location: " << point->requested_location
1333 out <<
"# DoF_index : Support location (for each component)\n";
1334 for (
unsigned int component = 0;
1335 component < dof_handler->get_fe(0).n_components();
1338 out << point->solution_indices[component] <<
" : "
1339 << point->support_point_locations[component] <<
'\n';
1346 out <<
"#Cannot access DoF_indices once cleared\n";
1351 if (independent_values.size() != 0)
1353 out <<
"Independent value(s): " << independent_values.size() <<
" : "
1354 << independent_values[0].size() <<
'\n';
1355 if (indep_names.size() > 0)
1358 for (
const auto &indep_name : indep_names)
1360 out <<
"<" << indep_name <<
"> ";
1367 out <<
"No independent values stored\n";
1370 if (data_store.begin() != data_store.end())
1373 <<
"Mnemonic: data set size (mask size, n true components) : n data sets\n";
1375 for (
const auto &data_entry : data_store)
1378 std::string vector_name = data_entry.first;
1379 typename std::map<std::string, ComponentMask>::iterator mask =
1380 component_mask.find(vector_name);
1381 Assert(mask != component_mask.end(),
1383 typename std::map<std::string, std::vector<std::string>>::iterator
1384 component_names = component_names_map.find(vector_name);
1385 Assert(component_names != component_names_map.end(),
1388 if (data_entry.second.size() != 0)
1390 out << data_entry.first <<
": " << data_entry.second.size() <<
" (";
1391 out << mask->second.size() <<
", "
1392 << mask->second.n_selected_components() <<
") : ";
1393 out << (data_entry.second)[0].size() <<
'\n';
1397 out << data_entry.first <<
": " << data_entry.second.size() <<
" (";
1398 out << mask->second.size() <<
", "
1399 << mask->second.n_selected_components() <<
") : ";
1400 out <<
"No points added" <<
'\n';
1403 if (component_names->second.size() > 0)
1405 for (
const auto &name : component_names->second)
1407 out <<
"<" << name <<
"> ";
1413 out <<
"***end of status output***\n\n";
1428 if (dataset_key.size() != independent_values[0].size())
1433 if (have_dof_handler)
1435 for (
const auto &data_entry : data_store)
1438 if ((data_entry.second)[0].size() != dataset_key.size())
1452 if (
std::abs(
static_cast<int>(dataset_key.size()) -
1453 static_cast<int>(independent_values[0].size())) >= 2)
1459 if (have_dof_handler)
1461 for (
const auto &data_entry : data_store)
1465 if (
std::abs(
static_cast<int>((data_entry.second)[0].size()) -
1466 static_cast<int>(dataset_key.size())) >= 2)
1495 triangulation_changed =
true;
1500#include "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
const std::vector< Point< spacedim > > & get_quadrature_points() const
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
const Point< spacedim > & quadrature_point(const unsigned int q) const
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
void add_field_name(const std::string &vector_name, const ComponentMask &component_mask=ComponentMask())
void evaluate_field_at_requested_location(const std::string &name, const VectorType &solution)
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)
void evaluate_field(const std::string &name, const VectorType &solution)
SmartPointer< const DoFHandler< dim >, PointValueHistory< dim > > dof_handler
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
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
@ 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.
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
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
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