40 namespace PointValueHistoryImplementation
46 const std::vector<types::global_dof_index> &new_sol_indices)
48 requested_location = new_requested_location;
49 support_point_locations = new_locations;
50 solution_indices = new_sol_indices;
59 const unsigned int n_independent_variables)
60 : n_indep(n_independent_variables)
72 std::vector<std::vector<double>>(
n_indep, std::vector<double>(0));
81 const unsigned int n_independent_variables)
82 : dof_handler(&dof_handler)
83 , n_indep(n_independent_variables)
95 std::vector<std::vector<double>>(
n_indep, std::vector<double>(0));
116 closed = point_value_history.
closed;
117 cleared = point_value_history.
cleared;
123 n_indep = point_value_history.
n_indep;
127 if (have_dof_handler)
130 dof_handler->get_triangulation().signals.any_change.connect(
131 [
this]() { this->tria_change_listener(); });
149 closed = point_value_history.
closed;
150 cleared = point_value_history.
cleared;
156 n_indep = point_value_history.
n_indep;
160 if (have_dof_handler)
163 dof_handler->get_triangulation().signals.any_change.connect(
164 [
this]() { this->tria_change_listener(); });
175 if (have_dof_handler)
177 tria_listener.disconnect();
191 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
192 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
204 dof_handler->get_fe().get_unit_support_points());
206 support_point_quadrature,
208 unsigned int n_support_points =
209 dof_handler->get_fe().get_unit_support_points().size();
210 unsigned int n_components = dof_handler->get_fe(0).n_components();
215 dof_handler->begin_active();
223 std::vector<unsigned int> current_fe_index(n_components,
226 std::vector<Point<dim>> current_points(n_components,
Point<dim>());
227 for (
unsigned int support_point = 0; support_point < n_support_points;
231 unsigned int component =
232 dof_handler->get_fe().system_to_component_index(support_point).first;
234 current_fe_index[component] = support_point;
247 for (; cell != endc; ++cell)
251 for (
unsigned int support_point = 0; support_point < n_support_points;
254 unsigned int component = dof_handler->get_fe()
255 .system_to_component_index(support_point)
261 location.
distance(current_points[component]))
264 current_points[component] = test_point;
266 current_fe_index[component] = support_point;
272 std::vector<types::global_dof_index> local_dof_indices(
273 dof_handler->get_fe().n_dofs_per_cell());
274 std::vector<types::global_dof_index> new_solution_indices;
275 current_cell->get_dof_indices(local_dof_indices);
295 for (
unsigned int component = 0;
296 component < dof_handler->get_fe(0).n_components();
299 new_solution_indices.push_back(
300 local_dof_indices[current_fe_index[component]]);
304 new_point_geometry_data(location, current_points, new_solution_indices);
305 point_geometry_data.push_back(new_point_geometry_data);
307 for (
auto &data_entry : data_store)
311 (component_mask.find(data_entry.first))->
second;
313 data_entry.second.resize(data_entry.second.size() + n_stored);
331 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
332 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
345 dof_handler->get_fe().get_unit_support_points());
347 support_point_quadrature,
349 unsigned int n_support_points =
350 dof_handler->get_fe().get_unit_support_points().size();
351 unsigned int n_components = dof_handler->get_fe(0).n_components();
356 dof_handler->begin_active();
366 std::vector<typename DoFHandler<dim>::active_cell_iterator> current_cell(
367 locations.size(), cell);
370 std::vector<Point<dim>> temp_points(n_components,
Point<dim>());
371 std::vector<unsigned int> temp_fe_index(n_components, 0);
372 for (
unsigned int support_point = 0; support_point < n_support_points;
376 unsigned int component =
377 dof_handler->get_fe().system_to_component_index(support_point).first;
379 temp_fe_index[component] = support_point;
381 std::vector<std::vector<Point<dim>>> current_points(
382 locations.size(), temp_points);
383 std::vector<std::vector<unsigned int>> current_fe_index(locations.size(),
395 for (; cell != endc; ++cell)
398 for (
unsigned int support_point = 0; support_point < n_support_points;
401 unsigned int component = dof_handler->get_fe()
402 .system_to_component_index(support_point)
407 for (
unsigned int point = 0; point < locations.size(); ++point)
409 if (locations[point].distance(test_point) <
410 locations[point].distance(current_points[point][component]))
413 current_points[point][component] = test_point;
414 current_cell[point] = cell;
415 current_fe_index[point][component] = support_point;
421 std::vector<types::global_dof_index> local_dof_indices(
422 dof_handler->get_fe().n_dofs_per_cell());
423 for (
unsigned int point = 0; point < locations.size(); ++point)
425 current_cell[point]->get_dof_indices(local_dof_indices);
426 std::vector<types::global_dof_index> new_solution_indices;
428 for (
unsigned int component = 0;
429 component < dof_handler->get_fe(0).n_components();
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);
441 point_geometry_data.push_back(new_point_geometry_data);
443 for (
auto &data_entry : data_store)
447 (component_mask.find(data_entry.first))->
second;
449 data_entry.second.resize(data_entry.second.size() + n_stored);
465 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
466 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
469 if (mask.represents_the_all_selected_mask() ==
false)
470 component_mask.insert(std::make_pair(vector_name, mask));
472 component_mask.insert(
473 std::make_pair(vector_name,
475 dof_handler->get_fe(0).n_components(),
true))));
480 std::pair<std::string, std::vector<std::string>> empty_names(
481 vector_name, std::vector<std::string>());
482 component_names_map.insert(empty_names);
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() :
491 dof_handler->get_fe(0).n_components());
494 point_geometry_data.size() * n_stored;
495 std::vector<std::vector<double>> vector_size(n_datastreams,
496 std::vector<double>(0));
497 pair_data.second = vector_size;
498 data_store.insert(pair_data);
505 const unsigned int n_components)
507 std::vector<bool> temp_mask(n_components,
true);
508 add_field_name(vector_name, temp_mask);
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 =
519 component_names_map.find(vector_name);
520 Assert(names != component_names_map.end(),
523 typename std::map<std::string, ComponentMask>::iterator mask =
524 component_mask.find(vector_name);
525 Assert(mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
526 unsigned int n_stored = mask->second.n_selected_components();
528 Assert(component_names.size() == n_stored,
531 names->second = component_names;
538 const std::vector<std::string> &independent_names)
540 Assert(independent_names.size() == n_indep,
543 indep_names = independent_names;
561 dof_handler =
nullptr;
562 have_dof_handler =
false;
578template <
typename VectorType>
581 const VectorType & solution)
587 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
588 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
594 static_cast<int>(independent_values[0].size())) < 2,
602 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
603 data_store_field = data_store.find(vector_name);
604 Assert(data_store_field != data_store.end(),
607 typename std::map<std::string, ComponentMask>::iterator mask =
608 component_mask.find(vector_name);
609 Assert(mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
611 unsigned int n_stored =
612 mask->second.n_selected_components(dof_handler->get_fe(0).n_components());
614 typename std::vector<
616 point = point_geometry_data.begin();
617 for (
unsigned int data_store_index = 0; point != point_geometry_data.end();
618 ++point, ++data_store_index)
625 for (
unsigned int store_index = 0, comp = 0;
626 comp < dof_handler->get_fe(0).n_components();
629 if (mask->second[comp])
631 unsigned int solution_index = point->solution_indices[comp];
633 ->second[data_store_index * n_stored + store_index]
646template <
typename VectorType>
649 const std::vector<std::string> &vector_names,
650 const VectorType & solution,
658 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
663 static_cast<int>(independent_values[0].size())) < 2,
673 "The update of normal vectors may not be requested for evaluation of "
674 "data on cells via DataPostprocessor."));
675 FEValues<dim> fe_values(dof_handler->get_fe(), quadrature, update_flags);
676 unsigned int n_components = dof_handler->get_fe(0).n_components();
677 unsigned int n_quadrature_points = quadrature.
size();
679 unsigned int n_output_variables = data_postprocessor.
get_names().size();
684 std::vector<typename VectorType::value_type> scalar_solution_values(
685 n_quadrature_points);
686 std::vector<Tensor<1, dim, typename VectorType::value_type>>
687 scalar_solution_gradients(n_quadrature_points);
688 std::vector<Tensor<2, dim, typename VectorType::value_type>>
689 scalar_solution_hessians(n_quadrature_points);
691 std::vector<Vector<typename VectorType::value_type>> vector_solution_values(
694 std::vector<std::vector<Tensor<1, dim, typename VectorType::value_type>>>
695 vector_solution_gradients(
700 std::vector<std::vector<Tensor<2, dim, typename VectorType::value_type>>>
701 vector_solution_hessians(
707 typename std::vector<
709 point = point_geometry_data.begin();
710 Assert(!dof_handler->get_triangulation().is_mixed_mesh(),
712 const auto reference_cell =
713 dof_handler->get_triangulation().get_reference_cells()[0];
714 for (
unsigned int data_store_index = 0; point != point_geometry_data.end();
715 ++point, ++data_store_index)
718 const Point<dim> requested_location = point->requested_location;
721 reference_cell.template get_default_linear_mapping<dim, dim>(),
728 std::vector<Vector<double>> computed_quantities(
732 std::vector<Point<dim>> quadrature_points =
734 double distance = cell->diameter();
735 unsigned int selected_point = 0;
736 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
738 if (requested_location.
distance(quadrature_points[q_point]) <
741 selected_point = q_point;
743 requested_location.
distance(quadrature_points[q_point]);
749 if (n_components == 1)
767 std::vector<double>(1, scalar_solution_values[selected_point]);
772 scalar_solution_gradients);
774 std::vector<Tensor<1, dim>>(
775 1, scalar_solution_gradients[selected_point]);
780 scalar_solution_hessians);
782 std::vector<Tensor<2, dim>>(
783 1, scalar_solution_hessians[selected_point]);
788 std::vector<Point<dim>>(1, quadrature_points[selected_point]);
792 computed_quantities);
804 std::copy(vector_solution_values[selected_point].begin(),
805 vector_solution_values[selected_point].end(),
811 vector_solution_gradients);
814 std::copy(vector_solution_gradients[selected_point].begin(),
815 vector_solution_gradients[selected_point].end(),
821 vector_solution_hessians);
824 std::copy(vector_solution_hessians[selected_point].begin(),
825 vector_solution_hessians[selected_point].end(),
830 std::vector<Point<dim>>(1, quadrature_points[selected_point]);
833 computed_quantities);
839 typename std::vector<std::string>::const_iterator name =
840 vector_names.begin();
841 for (; name != vector_names.end(); ++name)
843 typename std::map<std::string,
844 std::vector<std::vector<double>>>::iterator
845 data_store_field = data_store.find(*name);
846 Assert(data_store_field != data_store.end(),
849 typename std::map<std::string, ComponentMask>::iterator mask =
850 component_mask.find(*name);
851 Assert(mask != component_mask.end(),
854 unsigned int n_stored =
855 mask->second.n_selected_components(n_output_variables);
859 for (
unsigned int store_index = 0, comp = 0;
860 comp < n_output_variables;
863 if (mask->second[comp])
866 ->second[data_store_index * n_stored + store_index]
867 .push_back(computed_quantities[0](comp));
878template <
typename VectorType>
881 const std::string & vector_name,
882 const VectorType & solution,
886 std::vector<std::string> vector_names;
887 vector_names.push_back(vector_name);
888 evaluate_field(vector_names, solution, data_postprocessor, quadrature);
894template <
typename VectorType>
897 const std::string &vector_name,
898 const VectorType & solution)
900 using number =
typename VectorType::value_type;
905 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
911 static_cast<int>(independent_values[0].size())) < 2,
919 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
920 data_store_field = data_store.find(vector_name);
921 Assert(data_store_field != data_store.end(),
924 typename std::map<std::string, ComponentMask>::iterator mask =
925 component_mask.find(vector_name);
926 Assert(mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
928 unsigned int n_stored =
929 mask->second.n_selected_components(dof_handler->get_fe(0).n_components());
931 typename std::vector<
933 point = point_geometry_data.begin();
935 for (
unsigned int data_store_index = 0; point != point_geometry_data.end();
936 ++point, ++data_store_index)
943 point->requested_location,
948 for (
unsigned int store_index = 0, comp = 0; comp < mask->second.size();
951 if (mask->second[comp])
954 ->second[data_store_index * n_stored + store_index]
955 .push_back(value(comp));
971 Assert(deep_check(
false), ExcDataLostSync());
973 dataset_key.push_back(key);
981 const std::vector<double> &indep_values)
987 Assert(indep_values.size() == n_indep,
989 Assert(n_indep != 0, ExcNoIndependent());
991 static_cast<int>(independent_values[0].size())) < 2,
994 for (
unsigned int component = 0; component < n_indep; ++component)
995 independent_values[component].push_back(indep_values[component]);
1003 const std::string & base_name,
1004 const std::vector<
Point<dim>> &postprocessor_locations)
1013 std::string filename = base_name +
"_indep.gpl";
1014 std::ofstream to_gnuplot(filename.c_str());
1016 to_gnuplot <<
"# Data independent of mesh location\n";
1019 to_gnuplot <<
"# <Key> ";
1021 if (indep_names.size() > 0)
1023 for (
const auto &indep_name : indep_names)
1025 to_gnuplot <<
"<" << indep_name <<
"> ";
1031 for (
unsigned int component = 0; component < n_indep; ++component)
1033 to_gnuplot <<
"<Indep_" << component <<
"> ";
1038 for (
unsigned int key = 0; key < dataset_key.size(); ++key)
1040 to_gnuplot << dataset_key[key];
1042 for (
unsigned int component = 0; component < n_indep; ++component)
1044 to_gnuplot <<
" " << independent_values[component][key];
1055 if (have_dof_handler)
1057 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1058 AssertThrow(postprocessor_locations.size() == 0 ||
1059 postprocessor_locations.size() ==
1060 point_geometry_data.size(),
1062 point_geometry_data.size()));
1077 typename std::vector<internal::PointValueHistoryImplementation::
1078 PointGeometryData<dim>>::iterator point =
1079 point_geometry_data.begin();
1080 for (
unsigned int data_store_index = 0;
1081 point != point_geometry_data.end();
1082 ++point, ++data_store_index)
1086 std::string filename = base_name +
"_" +
1093 std::ofstream to_gnuplot(filename.c_str());
1098 to_gnuplot <<
"# Requested location: " << point->requested_location
1100 to_gnuplot <<
"# DoF_index : Support location (for each component)\n";
1101 for (
unsigned int component = 0;
1102 component < dof_handler->get_fe(0).n_components();
1105 to_gnuplot <<
"# " << point->solution_indices[component] <<
" : "
1106 << point->support_point_locations[component] <<
'\n';
1108 if (triangulation_changed)
1110 <<
"# (Original components and locations, may be invalidated by mesh change.)\n";
1112 if (postprocessor_locations.size() != 0)
1114 to_gnuplot <<
"# Postprocessor location: "
1115 << postprocessor_locations[data_store_index];
1116 if (triangulation_changed)
1117 to_gnuplot <<
" (may be approximate)\n";
1119 to_gnuplot <<
"#\n";
1123 to_gnuplot <<
"# <Key> ";
1125 if (indep_names.size() > 0)
1127 for (
const auto &indep_name : indep_names)
1129 to_gnuplot <<
"<" << indep_name <<
"> ";
1134 for (
unsigned int component = 0; component < n_indep; ++component)
1136 to_gnuplot <<
"<Indep_" << component <<
"> ";
1140 for (
const auto &data_entry : data_store)
1142 typename std::map<std::string, ComponentMask>::iterator mask =
1143 component_mask.find(data_entry.first);
1144 unsigned int n_stored = mask->second.n_selected_components();
1145 std::vector<std::string> names =
1146 (component_names_map.find(data_entry.first))->
second;
1148 if (names.size() > 0)
1152 for (
const auto &name : names)
1154 to_gnuplot <<
"<" << name <<
"> ";
1159 for (
unsigned int component = 0; component < n_stored;
1162 to_gnuplot <<
"<" << data_entry.first <<
"_" << component
1170 for (
unsigned int key = 0; key < dataset_key.size(); ++key)
1172 to_gnuplot << dataset_key[key];
1174 for (
unsigned int component = 0; component < n_indep; ++component)
1176 to_gnuplot <<
" " << independent_values[component][key];
1179 for (
const auto &data_entry : data_store)
1181 typename std::map<std::string, ComponentMask>::iterator mask =
1182 component_mask.find(data_entry.first);
1183 unsigned int n_stored = mask->second.n_selected_components();
1185 for (
unsigned int component = 0; component < n_stored;
1190 << (data_entry.second)[data_store_index * n_stored +
1211 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1212 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
1216 typename std::vector<
1218 point = point_geometry_data.begin();
1219 for (; point != point_geometry_data.end(); ++point)
1221 for (
unsigned int component = 0;
1222 component < dof_handler->get_fe(0).n_components();
1225 dof_vector(point->solution_indices[component]) = 1;
1235 std::vector<std::vector<
Point<dim>>> &locations)
1238 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1239 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
1241 std::vector<std::vector<Point<dim>>> actual_points;
1242 typename std::vector<
1244 point = point_geometry_data.begin();
1246 for (; point != point_geometry_data.end(); ++point)
1248 actual_points.push_back(point->support_point_locations);
1250 locations = actual_points;
1262 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1264 locations = std::vector<Point<dim>>();
1269 unsigned int n_quadrature_points = quadrature.
size();
1270 std::vector<Point<dim>> evaluation_points;
1273 typename std::vector<
1275 point = point_geometry_data.begin();
1276 Assert(!dof_handler->get_triangulation().is_mixed_mesh(),
1278 const auto reference_cell =
1279 dof_handler->get_triangulation().get_reference_cells()[0];
1280 for (
unsigned int data_store_index = 0; point != point_geometry_data.end();
1281 ++point, ++data_store_index)
1285 Point<dim> requested_location = point->requested_location;
1288 reference_cell.template get_default_linear_mapping<dim, dim>(),
1292 fe_values.reinit(cell);
1294 evaluation_points = fe_values.get_quadrature_points();
1295 double distance = cell->diameter();
1296 unsigned int selected_point = 0;
1298 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
1300 if (requested_location.
distance(evaluation_points[q_point]) <
1303 selected_point = q_point;
1305 requested_location.
distance(evaluation_points[q_point]);
1309 locations.push_back(evaluation_points[selected_point]);
1318 out <<
"***PointValueHistory status output***\n\n";
1319 out <<
"Closed: " << closed <<
'\n';
1320 out <<
"Cleared: " << cleared <<
'\n';
1321 out <<
"Triangulation_changed: " << triangulation_changed <<
'\n';
1322 out <<
"Have_dof_handler: " << have_dof_handler <<
'\n';
1323 out <<
"Geometric Data" <<
'\n';
1325 typename std::vector<
1327 point = point_geometry_data.begin();
1328 if (point == point_geometry_data.end())
1330 out <<
"No points stored currently\n";
1336 for (; point != point_geometry_data.end(); ++point)
1338 out <<
"# Requested location: " << point->requested_location
1340 out <<
"# DoF_index : Support location (for each component)\n";
1341 for (
unsigned int component = 0;
1342 component < dof_handler->get_fe(0).n_components();
1345 out << point->solution_indices[component] <<
" : "
1346 << point->support_point_locations[component] <<
'\n';
1353 out <<
"#Cannot access DoF_indices once cleared\n";
1358 if (independent_values.size() != 0)
1360 out <<
"Independent value(s): " << independent_values.size() <<
" : "
1361 << independent_values[0].size() <<
'\n';
1362 if (indep_names.size() > 0)
1365 for (
const auto &indep_name : indep_names)
1367 out <<
"<" << indep_name <<
"> ";
1374 out <<
"No independent values stored\n";
1377 if (data_store.begin() != data_store.end())
1380 <<
"Mnemonic: data set size (mask size, n true components) : n data sets\n";
1382 for (
const auto &data_entry : data_store)
1385 std::string vector_name = data_entry.first;
1386 typename std::map<std::string, ComponentMask>::iterator mask =
1387 component_mask.find(vector_name);
1388 Assert(mask != component_mask.end(),
1390 typename std::map<std::string, std::vector<std::string>>::iterator
1391 component_names = component_names_map.find(vector_name);
1392 Assert(component_names != component_names_map.end(),
1395 if (data_entry.second.size() != 0)
1397 out << data_entry.first <<
": " << data_entry.second.size() <<
" (";
1398 out << mask->second.size() <<
", "
1399 << mask->second.n_selected_components() <<
") : ";
1400 out << (data_entry.second)[0].size() <<
'\n';
1404 out << data_entry.first <<
": " << data_entry.second.size() <<
" (";
1405 out << mask->second.size() <<
", "
1406 << mask->second.n_selected_components() <<
") : ";
1407 out <<
"No points added" <<
'\n';
1410 if (component_names->second.size() > 0)
1412 for (
const auto &name : component_names->second)
1414 out <<
"<" << name <<
"> ";
1420 out <<
"***end of status output***\n\n";
1435 if (dataset_key.size() != independent_values[0].size())
1440 if (have_dof_handler)
1442 for (
const auto &data_entry : data_store)
1445 if ((data_entry.second)[0].size() != dataset_key.size())
1459 if (
std::abs(
static_cast<int>(dataset_key.size()) -
1460 static_cast<int>(independent_values[0].size())) >= 2)
1466 if (have_dof_handler)
1468 for (
const auto &data_entry : data_store)
1472 if (
std::abs(
static_cast<int>((data_entry.second)[0].size()) -
1473 static_cast<int>(dataset_key.size())) >= 2)
1502 triangulation_changed =
true;
1507#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
const Point< spacedim > & quadrature_point(const unsigned int q_point) const
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) 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
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
@ 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.
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