39 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,
225 fe_values.reinit(cell);
226 std::vector<Point<dim>> current_points(n_components,
Point<dim>());
227 for (
unsigned int support_point = 0; support_point < n_support_points;
232 unsigned int component =
233 dof_handler->get_fe().system_to_component_index(support_point).first;
234 current_points[component] = fe_values.quadrature_point(support_point);
235 current_fe_index[component] = support_point;
248 for (; cell != endc; cell++)
250 fe_values.reinit(cell);
252 for (
unsigned int support_point = 0; support_point < n_support_points;
255 unsigned int component = dof_handler->get_fe()
256 .system_to_component_index(support_point)
259 fe_values.quadrature_point(support_point);
262 location.
distance(current_points[component]))
265 current_points[component] = test_point;
267 current_fe_index[component] = support_point;
273 std::vector<types::global_dof_index> local_dof_indices(
274 dof_handler->get_fe().n_dofs_per_cell());
275 std::vector<types::global_dof_index> new_solution_indices;
276 current_cell->get_dof_indices(local_dof_indices);
296 for (
unsigned int component = 0;
297 component < dof_handler->get_fe(0).n_components();
300 new_solution_indices.push_back(
301 local_dof_indices[current_fe_index[component]]);
305 new_point_geometry_data(location, current_points, new_solution_indices);
306 point_geometry_data.push_back(new_point_geometry_data);
308 for (
auto &data_entry : data_store)
312 (component_mask.find(data_entry.first))->
second;
314 data_entry.second.resize(data_entry.second.size() + n_stored);
332 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
333 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
346 dof_handler->get_fe().get_unit_support_points());
348 support_point_quadrature,
350 unsigned int n_support_points =
351 dof_handler->get_fe().get_unit_support_points().size();
352 unsigned int n_components = dof_handler->get_fe(0).n_components();
357 dof_handler->begin_active();
367 std::vector<typename DoFHandler<dim>::active_cell_iterator> current_cell(
368 locations.size(), cell);
370 fe_values.reinit(cell);
371 std::vector<Point<dim>> temp_points(n_components,
Point<dim>());
372 std::vector<unsigned int> temp_fe_index(n_components, 0);
373 for (
unsigned int support_point = 0; support_point < n_support_points;
378 unsigned int component =
379 dof_handler->get_fe().system_to_component_index(support_point).first;
380 temp_points[component] = fe_values.quadrature_point(support_point);
381 temp_fe_index[component] = support_point;
383 std::vector<std::vector<Point<dim>>> current_points(
384 locations.size(), temp_points);
385 std::vector<std::vector<unsigned int>> current_fe_index(locations.size(),
397 for (; cell != endc; cell++)
399 fe_values.reinit(cell);
400 for (
unsigned int support_point = 0; support_point < n_support_points;
403 unsigned int component = dof_handler->get_fe()
404 .system_to_component_index(support_point)
407 fe_values.quadrature_point(support_point);
411 if (locations[
point].distance(test_point) <
412 locations[
point].distance(current_points[
point][component]))
415 current_points[
point][component] = test_point;
416 current_cell[
point] = cell;
417 current_fe_index[
point][component] = support_point;
423 std::vector<types::global_dof_index> local_dof_indices(
424 dof_handler->get_fe().n_dofs_per_cell());
427 current_cell[
point]->get_dof_indices(local_dof_indices);
428 std::vector<types::global_dof_index> new_solution_indices;
430 for (
unsigned int component = 0;
431 component < dof_handler->get_fe(0).n_components();
434 new_solution_indices.push_back(
435 local_dof_indices[current_fe_index[
point][component]]);
439 new_point_geometry_data(locations[
point],
440 current_points[
point],
441 new_solution_indices);
443 point_geometry_data.push_back(new_point_geometry_data);
445 for (
auto &data_entry : data_store)
449 (component_mask.find(data_entry.first))->
second;
451 data_entry.second.resize(data_entry.second.size() + n_stored);
467 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
468 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
472 component_mask.insert(std::make_pair(vector_name, mask));
474 component_mask.insert(
475 std::make_pair(vector_name,
477 dof_handler->get_fe(0).n_components(),
true))));
482 std::pair<std::string, std::vector<std::string>> empty_names(
483 vector_name, std::vector<std::string>());
484 component_names_map.insert(empty_names);
488 std::pair<std::string, std::vector<std::vector<double>>> pair_data;
489 pair_data.first = vector_name;
490 const unsigned int n_stored =
493 dof_handler->get_fe(0).n_components());
496 point_geometry_data.size() * n_stored;
497 std::vector<std::vector<double>> vector_size(n_datastreams,
498 std::vector<double>(0));
499 pair_data.second = vector_size;
500 data_store.insert(pair_data);
507 const unsigned int n_components)
509 std::vector<bool> temp_mask(n_components,
true);
510 add_field_name(vector_name, temp_mask);
517 const std::string & vector_name,
518 const std::vector<std::string> &component_names)
520 typename std::map<std::string, std::vector<std::string>>::iterator names =
521 component_names_map.find(vector_name);
522 Assert(names != component_names_map.end(),
525 typename std::map<std::string, ComponentMask>::iterator mask =
526 component_mask.find(vector_name);
527 Assert(mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
528 unsigned int n_stored = mask->second.n_selected_components();
530 Assert(component_names.size() == n_stored,
533 names->second = component_names;
540 const std::vector<std::string> &independent_names)
542 Assert(independent_names.size() == n_indep,
545 indep_names = independent_names;
563 dof_handler =
nullptr;
564 have_dof_handler =
false;
580template <
typename VectorType>
583 const VectorType & solution)
589 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
590 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
596 static_cast<int>(independent_values[0].size())) < 2,
604 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
605 data_store_field = data_store.find(vector_name);
606 Assert(data_store_field != data_store.end(),
609 typename std::map<std::string, ComponentMask>::iterator mask =
610 component_mask.find(vector_name);
611 Assert(mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
613 unsigned int n_stored =
614 mask->second.n_selected_components(dof_handler->get_fe(0).n_components());
616 typename std::vector<
618 point = point_geometry_data.begin();
619 for (
unsigned int data_store_index = 0;
point != point_geometry_data.end();
620 ++
point, ++data_store_index)
627 for (
unsigned int store_index = 0, comp = 0;
628 comp < dof_handler->get_fe(0).n_components();
631 if (mask->second[comp])
633 unsigned int solution_index =
point->solution_indices[comp];
635 ->second[data_store_index * n_stored + store_index]
648template <
typename VectorType>
651 const std::vector<std::string> &vector_names,
652 const VectorType & solution,
660 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
665 static_cast<int>(independent_values[0].size())) < 2,
675 "The update of normal vectors may not be requested for evaluation of "
676 "data on cells via DataPostprocessor."));
677 FEValues<dim> fe_values(dof_handler->get_fe(), quadrature, update_flags);
678 unsigned int n_components = dof_handler->get_fe(0).n_components();
679 unsigned int n_quadrature_points = quadrature.
size();
681 unsigned int n_output_variables = data_postprocessor.
get_names().size();
686 std::vector<typename VectorType::value_type> scalar_solution_values(
687 n_quadrature_points);
688 std::vector<Tensor<1, dim, typename VectorType::value_type>>
689 scalar_solution_gradients(n_quadrature_points);
690 std::vector<Tensor<2, dim, typename VectorType::value_type>>
691 scalar_solution_hessians(n_quadrature_points);
693 std::vector<Vector<typename VectorType::value_type>> vector_solution_values(
696 std::vector<std::vector<Tensor<1, dim, typename VectorType::value_type>>>
697 vector_solution_gradients(
702 std::vector<std::vector<Tensor<2, dim, typename VectorType::value_type>>>
703 vector_solution_hessians(
709 typename std::vector<
711 point = point_geometry_data.begin();
712 for (
unsigned int data_store_index = 0;
point != point_geometry_data.end();
713 ++
point, ++data_store_index)
724 fe_values.reinit(cell);
725 std::vector<Vector<double>> computed_quantities(
730 fe_values.get_quadrature_points();
731 double distance = cell->diameter();
732 unsigned int selected_point = 0;
733 for (
unsigned int q_point = 0; q_point < n_quadrature_points; q_point++)
738 selected_point = q_point;
746 if (n_components == 1)
762 fe_values.get_function_values(solution, scalar_solution_values);
764 std::vector<double>(1, scalar_solution_values[selected_point]);
768 fe_values.get_function_gradients(solution,
769 scalar_solution_gradients);
771 std::vector<Tensor<1, dim>>(
772 1, scalar_solution_gradients[selected_point]);
776 fe_values.get_function_hessians(solution,
777 scalar_solution_hessians);
779 std::vector<Tensor<2, dim>>(
780 1, scalar_solution_hessians[selected_point]);
789 computed_quantities);
798 fe_values.get_function_values(solution, vector_solution_values);
802 vector_solution_values[selected_point].
end(),
807 fe_values.get_function_gradients(solution,
808 vector_solution_gradients);
812 vector_solution_gradients[selected_point].
end(),
817 fe_values.get_function_hessians(solution,
818 vector_solution_hessians);
822 vector_solution_hessians[selected_point].
end(),
830 computed_quantities);
836 typename std::vector<std::string>::const_iterator name =
837 vector_names.begin();
838 for (; name != vector_names.end(); ++name)
840 typename std::map<std::string,
841 std::vector<std::vector<double>>>::iterator
842 data_store_field = data_store.find(*name);
843 Assert(data_store_field != data_store.end(),
846 typename std::map<std::string, ComponentMask>::iterator mask =
847 component_mask.find(*name);
848 Assert(mask != component_mask.end(),
851 unsigned int n_stored =
852 mask->second.n_selected_components(n_output_variables);
856 for (
unsigned int store_index = 0, comp = 0;
857 comp < n_output_variables;
860 if (mask->second[comp])
863 ->second[data_store_index * n_stored + store_index]
864 .push_back(computed_quantities[0](comp));
875template <
typename VectorType>
878 const std::string & vector_name,
879 const VectorType & solution,
883 std::vector<std::string> vector_names;
884 vector_names.push_back(vector_name);
885 evaluate_field(vector_names, solution, data_postprocessor, quadrature);
891template <
typename VectorType>
894 const std::string &vector_name,
895 const VectorType & solution)
897 using number =
typename VectorType::value_type;
902 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
908 static_cast<int>(independent_values[0].size())) < 2,
916 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
917 data_store_field = data_store.find(vector_name);
918 Assert(data_store_field != data_store.end(),
921 typename std::map<std::string, ComponentMask>::iterator mask =
922 component_mask.find(vector_name);
923 Assert(mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
925 unsigned int n_stored =
926 mask->second.n_selected_components(dof_handler->get_fe(0).n_components());
928 typename std::vector<
930 point = point_geometry_data.begin();
932 for (
unsigned int data_store_index = 0;
point != point_geometry_data.end();
933 ++
point, ++data_store_index)
940 point->requested_location,
945 for (
unsigned int store_index = 0, comp = 0; comp < mask->second.size();
948 if (mask->second[comp])
951 ->second[data_store_index * n_stored + store_index]
952 .push_back(value(comp));
968 Assert(deep_check(
false), ExcDataLostSync());
970 dataset_key.push_back(key);
978 const std::vector<double> &indep_values)
984 Assert(indep_values.size() == n_indep,
986 Assert(n_indep != 0, ExcNoIndependent());
988 static_cast<int>(independent_values[0].size())) < 2,
991 for (
unsigned int component = 0; component < n_indep; component++)
992 independent_values[component].push_back(indep_values[component]);
1000 const std::string & base_name,
1001 const std::vector<
Point<dim>> &postprocessor_locations)
1010 std::string filename = base_name +
"_indep.gpl";
1011 std::ofstream to_gnuplot(filename.c_str());
1013 to_gnuplot <<
"# Data independent of mesh location\n";
1016 to_gnuplot <<
"# <Key> ";
1018 if (indep_names.size() > 0)
1020 for (
const auto &indep_name : indep_names)
1022 to_gnuplot <<
"<" << indep_name <<
"> ";
1028 for (
unsigned int component = 0; component < n_indep; component++)
1030 to_gnuplot <<
"<Indep_" << component <<
"> ";
1035 for (
unsigned int key = 0; key < dataset_key.size(); key++)
1037 to_gnuplot << dataset_key[key];
1039 for (
unsigned int component = 0; component < n_indep; component++)
1041 to_gnuplot <<
" " << independent_values[component][key];
1052 if (have_dof_handler)
1054 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1055 AssertThrow(postprocessor_locations.size() == 0 ||
1056 postprocessor_locations.size() ==
1057 point_geometry_data.size(),
1059 point_geometry_data.size()));
1074 typename std::vector<internal::PointValueHistoryImplementation::
1075 PointGeometryData<dim>>::iterator
point =
1076 point_geometry_data.begin();
1077 for (
unsigned int data_store_index = 0;
1078 point != point_geometry_data.end();
1079 ++
point, ++data_store_index)
1083 std::string filename = base_name +
"_" +
1090 std::ofstream to_gnuplot(filename.c_str());
1095 to_gnuplot <<
"# Requested location: " <<
point->requested_location
1097 to_gnuplot <<
"# DoF_index : Support location (for each component)\n";
1098 for (
unsigned int component = 0;
1099 component < dof_handler->get_fe(0).n_components();
1102 to_gnuplot <<
"# " <<
point->solution_indices[component] <<
" : "
1103 <<
point->support_point_locations[component] <<
"\n";
1105 if (triangulation_changed)
1107 <<
"# (Original components and locations, may be invalidated by mesh change.)\n";
1109 if (postprocessor_locations.size() != 0)
1111 to_gnuplot <<
"# Postprocessor location: "
1112 << postprocessor_locations[data_store_index];
1113 if (triangulation_changed)
1114 to_gnuplot <<
" (may be approximate)\n";
1116 to_gnuplot <<
"#\n";
1120 to_gnuplot <<
"# <Key> ";
1122 if (indep_names.size() > 0)
1124 for (
const auto &indep_name : indep_names)
1126 to_gnuplot <<
"<" << indep_name <<
"> ";
1131 for (
unsigned int component = 0; component < n_indep; component++)
1133 to_gnuplot <<
"<Indep_" << component <<
"> ";
1137 for (
const auto &data_entry : data_store)
1139 typename std::map<std::string, ComponentMask>::iterator mask =
1140 component_mask.find(data_entry.first);
1141 unsigned int n_stored = mask->second.n_selected_components();
1142 std::vector<std::string> names =
1143 (component_names_map.find(data_entry.first))->
second;
1145 if (names.size() > 0)
1149 for (
const auto &name : names)
1151 to_gnuplot <<
"<" << name <<
"> ";
1156 for (
unsigned int component = 0; component < n_stored;
1159 to_gnuplot <<
"<" << data_entry.first <<
"_" << component
1167 for (
unsigned int key = 0; key < dataset_key.size(); key++)
1169 to_gnuplot << dataset_key[key];
1171 for (
unsigned int component = 0; component < n_indep; component++)
1173 to_gnuplot <<
" " << independent_values[component][key];
1176 for (
const auto &data_entry : data_store)
1178 typename std::map<std::string, ComponentMask>::iterator mask =
1179 component_mask.find(data_entry.first);
1180 unsigned int n_stored = mask->second.n_selected_components();
1182 for (
unsigned int component = 0; component < n_stored;
1187 << (data_entry.second)[data_store_index * n_stored +
1208 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1209 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
1213 typename std::vector<
1215 point = point_geometry_data.begin();
1216 for (;
point != point_geometry_data.end(); ++
point)
1218 for (
unsigned int component = 0;
1219 component < dof_handler->get_fe(0).n_components();
1222 dof_vector(
point->solution_indices[component]) = 1;
1232 std::vector<std::vector<
Point<dim>>> &locations)
1235 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1236 AssertThrow(!triangulation_changed, ExcDoFHandlerChanged());
1238 std::vector<std::vector<Point<dim>>> actual_points;
1239 typename std::vector<
1241 point = point_geometry_data.begin();
1243 for (;
point != point_geometry_data.end(); ++
point)
1245 actual_points.push_back(
point->support_point_locations);
1247 locations = actual_points;
1259 AssertThrow(have_dof_handler, ExcDoFHandlerRequired());
1261 locations = std::vector<Point<dim>>();
1266 unsigned int n_quadrature_points = quadrature.
size();
1267 std::vector<Point<dim>> evaluation_points;
1270 typename std::vector<
1272 point = point_geometry_data.begin();
1273 for (
unsigned int data_store_index = 0;
point != point_geometry_data.end();
1274 ++
point, ++data_store_index)
1284 fe_values.reinit(cell);
1286 evaluation_points = fe_values.get_quadrature_points();
1287 double distance = cell->diameter();
1288 unsigned int selected_point = 0;
1290 for (
unsigned int q_point = 0; q_point < n_quadrature_points; q_point++)
1292 if (requested_location.
distance(evaluation_points[q_point]) <
1295 selected_point = q_point;
1297 requested_location.
distance(evaluation_points[q_point]);
1301 locations.push_back(evaluation_points[selected_point]);
1310 out <<
"***PointValueHistory status output***\n\n";
1311 out <<
"Closed: " << closed <<
"\n";
1312 out <<
"Cleared: " << cleared <<
"\n";
1313 out <<
"Triangulation_changed: " << triangulation_changed <<
"\n";
1314 out <<
"Have_dof_handler: " << have_dof_handler <<
"\n";
1315 out <<
"Geometric Data"
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"
1404 if (component_names->second.size() > 0)
1406 for (
const auto &name : component_names->second)
1408 out <<
"<" << name <<
"> ";
1414 out <<
"***end of status output***\n\n";
1429 if (dataset_key.size() != independent_values[0].size())
1434 if (have_dof_handler)
1436 for (
const auto &data_entry : data_store)
1439 if ((data_entry.second)[0].size() != dataset_key.size())
1453 if (
std::abs(
static_cast<int>(dataset_key.size()) -
1454 static_cast<int>(independent_values[0].size())) >= 2)
1460 if (have_dof_handler)
1462 for (
const auto &data_entry : data_store)
1466 if (
std::abs(
static_cast<int>((data_entry.second)[0].size()) -
1467 static_cast<int>(dataset_key.size())) >= 2)
1496 triangulation_changed =
true;
1501#include "point_value_history.inst"
bool represents_the_all_selected_mask() const
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 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)
Only a constructor needed for this class (a struct really)
#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
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
void quadrature_points(const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const std::vector< std::vector< double > > &properties={})
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
void copy(const T *begin, const T *end, U *dest)
::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