17 #include <deal.II/lac/vector.h> 18 #include <deal.II/lac/block_vector.h> 19 #include <deal.II/lac/la_vector.h> 20 #include <deal.II/lac/la_parallel_vector.h> 21 #include <deal.II/lac/la_parallel_block_vector.h> 22 #include <deal.II/lac/petsc_parallel_vector.h> 23 #include <deal.II/lac/petsc_parallel_block_vector.h> 24 #include <deal.II/lac/trilinos_vector.h> 25 #include <deal.II/lac/trilinos_parallel_block_vector.h> 26 #include <deal.II/lac/vector_element_access.h> 28 #include <deal.II/numerics/vector_tools.h> 30 #include <deal.II/numerics/point_value_history.h> 35 DEAL_II_NAMESPACE_OPEN
40 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;
61 n_indep (n_independent_variables)
65 triangulation_changed =
false;
66 have_dof_handler =
false;
69 dataset_key = std::vector <double> ();
73 = std::vector<std::vector <double> > (n_indep, std::vector <double> (0));
74 indep_names = std::vector <std::string> ();
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));
115 closed = point_value_history.
closed;
116 cleared = point_value_history.
cleared;
122 n_indep = point_value_history.
n_indep;
126 if (have_dof_handler)
147 closed = point_value_history.
closed;
148 cleared = point_value_history.
cleared;
154 n_indep = point_value_history.
n_indep;
158 if (have_dof_handler)
173 if (have_dof_handler)
175 tria_listener.disconnect ();
189 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
190 AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
202 support_point_quadrature (dof_handler->get_fe().get_unit_support_points ());
204 support_point_quadrature,
206 unsigned int n_support_points
207 = dof_handler->
get_fe().get_unit_support_points ().size ();
209 = dof_handler->get_fe(0).n_components ();
216 endc = dof_handler->
end ();
223 std::vector <unsigned int> current_fe_index (
n_components, 0);
224 fe_values.reinit (cell);
226 for (
unsigned int support_point = 0;
227 support_point < n_support_points; support_point++)
231 unsigned int component
232 = dof_handler->get_fe().system_to_component_index (support_point).first;
233 current_points [component] = fe_values.quadrature_point (support_point);
234 current_fe_index [component] = support_point;
247 for (; cell != endc; cell++)
249 fe_values.reinit (cell);
251 for (
unsigned int support_point = 0;
252 support_point < n_support_points; support_point++)
254 unsigned int component
255 = dof_handler->get_fe().system_to_component_index (support_point).first;
257 = fe_values.quadrature_point (support_point);
259 if (location.
distance (test_point) <
260 location.
distance (current_points [component]))
263 current_points [component] = test_point;
265 current_fe_index [component] = support_point;
271 std::vector<types::global_dof_index>
272 local_dof_indices (dof_handler->get_fe().dofs_per_cell);
273 std::vector <types::global_dof_index> new_solution_indices;
274 current_cell->get_dof_indices (local_dof_indices);
294 for (
unsigned int component = 0;
295 component < dof_handler->get_fe(0).n_components (); component++)
298 .push_back (local_dof_indices[current_fe_index [component]]);
302 new_point_geometry_data (location, current_points, new_solution_indices);
303 point_geometry_data.push_back (new_point_geometry_data);
305 std::map <std::string, std::vector <std::vector <double> > >::iterator
306 data_store_begin = data_store.begin ();
307 for (; data_store_begin != data_store.end (); ++data_store_begin)
311 const ComponentMask ¤t_mask = (component_mask.find (data_store_begin->first))->second;
313 data_store_begin->second.resize(data_store_begin->second.size() + n_stored);
331 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
332 AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
344 Quadrature<dim> support_point_quadrature (dof_handler->get_fe().get_unit_support_points ());
346 unsigned int n_support_points = dof_handler->
get_fe().get_unit_support_points ().size ();
347 unsigned int n_components = dof_handler->get_fe(0).n_components ();
361 std::vector <typename DoFHandler<dim>::active_cell_iterator > current_cell (locations.size (), cell);
363 fe_values.reinit (cell);
365 std::vector <unsigned int> temp_fe_index (
n_components, 0);
366 for (
unsigned int support_point = 0; support_point < n_support_points; support_point++)
370 unsigned int component = dof_handler->get_fe().system_to_component_index (support_point).first;
371 temp_points [component] = fe_values.quadrature_point (support_point);
372 temp_fe_index [component] = support_point;
374 std::vector <std::vector <Point <dim> > > current_points (locations.size (), temp_points);
375 std::vector <std::vector <unsigned int> > current_fe_index (locations.size (), temp_fe_index);
386 for (; cell != endc; cell++)
388 fe_values.reinit (cell);
389 for (
unsigned int support_point = 0; support_point < n_support_points; support_point++)
391 unsigned int component = dof_handler->get_fe().system_to_component_index (support_point).first;
392 Point<dim> test_point = fe_values.quadrature_point (support_point);
394 for (
unsigned int point = 0; point < locations.size (); point++)
396 if (locations[point].distance (test_point) < locations[point].distance (current_points[point][component]))
399 current_points[point][component] = test_point;
400 current_cell[point] = cell;
401 current_fe_index[point][component] = support_point;
407 std::vector<types::global_dof_index> local_dof_indices (dof_handler->get_fe().dofs_per_cell);
408 for (
unsigned int point = 0; point < locations.size (); point++)
410 current_cell[point]->get_dof_indices (local_dof_indices);
411 std::vector<types::global_dof_index> new_solution_indices;
413 for (
unsigned int component = 0; component < dof_handler->get_fe(0).n_components (); component++)
415 new_solution_indices.push_back (local_dof_indices[current_fe_index[point][component]]);
420 point_geometry_data.push_back (new_point_geometry_data);
422 std::map <std::string, std::vector <std::vector <double> > >::iterator
423 data_store_begin = data_store.begin ();
424 for (; data_store_begin != data_store.end (); ++data_store_begin)
428 const ComponentMask current_mask = (component_mask.find (data_store_begin->first))->second;
430 data_store_begin->second.resize(data_store_begin->second.size() + n_stored);
449 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
450 AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
454 component_mask.insert (std::make_pair (vector_name, mask));
456 component_mask.insert (std::make_pair (vector_name,
457 ComponentMask(std::vector<bool>(dof_handler->get_fe(0).n_components(),
true))));
462 std::pair <std::string, std::vector <std::string> >
463 empty_names (vector_name, std::vector <std::string> ());
464 component_names_map.insert (empty_names);
468 std::pair<std::string, std::vector <std::vector <double> > > pair_data;
469 pair_data.first = vector_name;
474 dof_handler->get_fe(0).n_components());
476 int n_datastreams = point_geometry_data.size () * n_stored;
477 std::vector < std::vector <double> > vector_size (n_datastreams,
478 std::vector <double> (0));
479 pair_data.second = vector_size;
480 data_store.insert (pair_data);
489 add_field_name (vector_name, temp_mask);
496 const std::vector <std::string> &component_names)
498 typename std::map <std::string, std::vector <std::string> >::iterator names = component_names_map.find(vector_name);
499 Assert (names != component_names_map.end(),
ExcMessage(
"vector_name not in class"));
501 typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(vector_name);
502 Assert (mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
503 unsigned int n_stored = mask->second.n_selected_components();
507 names->second = component_names;
517 indep_names = independent_names;
535 dof_handler =
nullptr;
536 have_dof_handler =
false;
552 template <
typename VectorType>
555 const VectorType &solution)
561 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
562 AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
566 Assert (std::abs ((
int) dataset_key.size () - (int) independent_values[0].size ()) < 2, ExcDataLostSync ());
573 typename std::map <std::string, std::vector <std::vector <double> > >::iterator data_store_field = data_store.find(vector_name);
574 Assert (data_store_field != data_store.end(),
ExcMessage(
"vector_name not in class"));
576 typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(vector_name);
577 Assert (mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
579 unsigned int n_stored = mask->second.n_selected_components(dof_handler->get_fe(0).n_components ());
581 typename std::vector <internal::PointValueHistoryImplementation::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
582 for (
unsigned int data_store_index = 0; point != point_geometry_data.end (); ++point, ++data_store_index)
589 for (
unsigned int store_index = 0, comp = 0; comp < dof_handler->get_fe(0).n_components (); comp++)
591 if (mask->second[comp])
593 unsigned int solution_index = point->solution_indices[comp];
594 data_store_field->second[data_store_index * n_stored + store_index].push_back (
595 internal::ElementAccess<VectorType>::get(solution,solution_index));
607 template <
typename VectorType>
610 const VectorType &solution,
618 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
621 Assert (std::abs ((
int) dataset_key.size () - (int) independent_values[0].size ()) < 2, ExcDataLostSync ());
627 ExcMessage(
"The update of normal vectors may not be requested for evaluation of " 628 "data on cells via DataPostprocessor."));
629 FEValues<dim> fe_values (dof_handler->get_fe(), quadrature, update_flags);
630 unsigned int n_components = dof_handler->get_fe(0).n_components ();
631 unsigned int n_quadrature_points = quadrature.
size();
633 unsigned int n_output_variables = data_postprocessor.
get_names().size();
638 std::vector<typename VectorType::value_type > scalar_solution_values (n_quadrature_points);
639 std::vector<Tensor< 1, dim, typename VectorType::value_type > > scalar_solution_gradients (n_quadrature_points);
640 std::vector<Tensor< 2, dim, typename VectorType::value_type > > scalar_solution_hessians (n_quadrature_points);
642 std::vector< Vector< typename VectorType::value_type > >
645 std::vector< std::vector< Tensor< 1, dim, typename VectorType::value_type > > >
646 vector_solution_gradients (n_quadrature_points,
650 std::vector< std::vector< Tensor< 2, dim, typename VectorType::value_type > > >
651 vector_solution_hessians (n_quadrature_points,
656 typename std::vector <internal::PointValueHistoryImplementation::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
657 for (
unsigned int data_store_index = 0; point != point_geometry_data.end (); ++point, ++data_store_index)
660 const Point <dim> requested_location = point->requested_location;
665 fe_values.reinit (cell);
666 std::vector< Vector< double > > computed_quantities (1,
Vector <double> (n_output_variables));
669 std::vector<Point<dim> > quadrature_points = fe_values.get_quadrature_points();
670 double distance = cell->diameter ();
671 unsigned int selected_point = 0;
672 for (
unsigned int q_point = 0; q_point < n_quadrature_points; q_point++)
674 if (requested_location.
distance (quadrature_points[q_point]) < distance)
676 selected_point = q_point;
677 distance = requested_location.
distance (quadrature_points[q_point]);
699 fe_values.get_function_values (solution,
700 scalar_solution_values);
702 = std::vector<double> (1, scalar_solution_values[selected_point]);
706 fe_values.get_function_gradients (solution,
707 scalar_solution_gradients);
709 = std::vector<Tensor<1,dim>> (1, scalar_solution_gradients[selected_point]);
713 fe_values.get_function_hessians (solution,
714 scalar_solution_hessians);
716 = std::vector<Tensor<2,dim>> (1, scalar_solution_hessians[selected_point]);
721 = std::vector<Point<dim>>(1, quadrature_points[selected_point]);
725 computed_quantities);
734 fe_values.get_function_values (solution,
735 vector_solution_values);
737 std::copy (vector_solution_values[selected_point].begin(),
738 vector_solution_values[selected_point].end(),
743 fe_values.get_function_gradients (solution,
744 vector_solution_gradients);
746 std::copy (vector_solution_gradients[selected_point].begin(),
747 vector_solution_gradients[selected_point].end(),
752 fe_values.get_function_hessians (solution,
753 vector_solution_hessians);
755 std::copy (vector_solution_hessians[selected_point].begin(),
756 vector_solution_hessians[selected_point].end(),
761 = std::vector<Point<dim>>(1, quadrature_points[selected_point]);
764 computed_quantities);
770 typename std::vector<std::string>::const_iterator name = vector_names.begin();
771 for (; name != vector_names.end(); ++name)
773 typename std::map <std::string, std::vector <std::vector <double> > >::iterator data_store_field = data_store.find(*name);
774 Assert (data_store_field != data_store.end(),
ExcMessage(
"vector_name not in class"));
776 typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(*name);
777 Assert (mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
779 unsigned int n_stored = mask->second.n_selected_components(n_output_variables);
783 for (
unsigned int store_index = 0, comp = 0; comp < n_output_variables; comp++)
785 if (mask->second[comp])
787 data_store_field->second[data_store_index * n_stored + store_index].push_back (computed_quantities[0](comp));
798 template <
typename VectorType>
801 const VectorType &solution,
805 std::vector <std::string> vector_names;
806 vector_names.push_back (vector_name);
807 evaluate_field (vector_names, solution, data_postprocessor, quadrature);
813 template <
typename VectorType>
816 const VectorType &solution)
818 typedef typename VectorType::value_type number;
823 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
827 Assert (std::abs ((
int) dataset_key.size () - (int) independent_values[0].size ()) < 2, ExcDataLostSync ());
834 typename std::map <std::string, std::vector <std::vector <double> > >::iterator data_store_field = data_store.find(vector_name);
835 Assert (data_store_field != data_store.end(),
ExcMessage(
"vector_name not in class"));
837 typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(vector_name);
838 Assert (mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
840 unsigned int n_stored = mask->second.n_selected_components(dof_handler->get_fe(0).n_components ());
842 typename std::vector <internal::PointValueHistoryImplementation::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
844 for (
unsigned int data_store_index = 0; point != point_geometry_data.end (); ++point, ++data_store_index)
853 for (
unsigned int store_index = 0, comp = 0; comp < mask->second.size(); comp++)
855 if (mask->second[comp])
857 data_store_field->second[data_store_index * n_stored + store_index].push_back (value (comp));
873 Assert (deep_check (
false), ExcDataLostSync ());
875 dataset_key.push_back (key);
889 Assert (n_indep != 0, ExcNoIndependent ());
890 Assert (std::abs ((
int) dataset_key.size () - (int) independent_values[0].size ()) < 2, ExcDataLostSync ());
892 for (
unsigned int component = 0; component < n_indep; component++)
893 independent_values[component].push_back (indep_values[component]);
901 const std::vector <
Point <dim> > &postprocessor_locations)
905 AssertThrow (deep_check (
true), ExcDataLostSync ());
910 std::string filename = base_name +
"_indep.gpl";
911 std::ofstream to_gnuplot (filename.c_str ());
913 to_gnuplot <<
"# Data independent of mesh location\n";
916 to_gnuplot <<
"# <Key> ";
918 if (indep_names.size() > 0)
920 for (
unsigned int name = 0; name < indep_names.size(); name++)
922 to_gnuplot <<
"<" << indep_names [name] <<
"> ";
928 for (
unsigned int component = 0; component < n_indep; component++)
930 to_gnuplot <<
"<Indep_" << component <<
"> ";
935 for (
unsigned int key = 0; key < dataset_key.size (); key++)
937 to_gnuplot << dataset_key[key];
939 for (
unsigned int component = 0; component < n_indep; component++)
941 to_gnuplot <<
" " << independent_values[component][key];
952 if (have_dof_handler)
954 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
955 AssertThrow (postprocessor_locations.size() == 0 || postprocessor_locations.size() == point_geometry_data.size(),
ExcDimensionMismatch (postprocessor_locations.size(), point_geometry_data.size()));
970 typename std::vector <internal::PointValueHistoryImplementation::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
971 for (
unsigned int data_store_index = 0; point != point_geometry_data.end (); ++point, ++data_store_index)
980 std::ofstream to_gnuplot (filename.c_str ());
985 to_gnuplot <<
"# Requested location: " << point->requested_location <<
"\n";
986 to_gnuplot <<
"# DoF_index : Support location (for each component)\n";
987 for (
unsigned int component = 0; component < dof_handler->get_fe(0).n_components (); component++)
989 to_gnuplot <<
"# " << point->solution_indices[component] <<
" : " << point->support_point_locations [component] <<
"\n";
991 if (triangulation_changed)
992 to_gnuplot <<
"# (Original components and locations, may be invalidated by mesh change.)\n";
994 if (postprocessor_locations.size() != 0)
996 to_gnuplot <<
"# Postprocessor location: " << postprocessor_locations[data_store_index];
997 if (triangulation_changed)
998 to_gnuplot <<
" (may be approximate)\n";
1000 to_gnuplot <<
"#\n";
1004 to_gnuplot <<
"# <Key> ";
1006 if (indep_names.size() > 0)
1008 for (
unsigned int name = 0; name < indep_names.size(); name++)
1010 to_gnuplot <<
"<" << indep_names [name] <<
"> ";
1015 for (
unsigned int component = 0; component < n_indep; component++)
1017 to_gnuplot <<
"<Indep_" << component <<
"> ";
1021 for (std::map <std::string, std::vector <std::vector <double> > >::iterator
1022 data_store_begin = data_store.begin (); data_store_begin != data_store.end (); ++data_store_begin)
1024 typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(data_store_begin->first);
1025 unsigned int n_stored = mask->second.n_selected_components();
1026 std::vector <std::string> names = (component_names_map.find (data_store_begin->first))->second;
1028 if (names.size() > 0)
1031 for (
unsigned int component = 0; component < names.size(); component++)
1033 to_gnuplot <<
"<" << names[component] <<
"> ";
1038 for (
unsigned int component = 0; component < n_stored; component++)
1040 to_gnuplot <<
"<" << data_store_begin->first <<
"_" << component <<
"> ";
1047 for (
unsigned int key = 0; key < dataset_key.size (); key++)
1049 to_gnuplot << dataset_key[key];
1051 for (
unsigned int component = 0; component < n_indep; component++)
1053 to_gnuplot <<
" " << independent_values[component][key];
1056 for (std::map <std::string, std::vector <std::vector <double> > >::iterator
1057 data_store_begin = data_store.begin ();
1058 data_store_begin != data_store.end (); ++data_store_begin)
1060 typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(data_store_begin->first);
1061 unsigned int n_stored = mask->second.n_selected_components();
1063 for (
unsigned int component = 0; component < n_stored; component++)
1065 to_gnuplot <<
" " << (data_store_begin->second)[data_store_index * n_stored + component][key];
1071 to_gnuplot.close ();
1085 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
1086 AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
1090 typename std::vector <internal::PointValueHistoryImplementation::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
1091 for (; point != point_geometry_data.end (); ++point)
1093 for (
unsigned int component = 0; component < dof_handler->get_fe(0).n_components (); component++)
1095 dof_vector (point->solution_indices[component]) = 1;
1107 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
1108 AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
1110 std::vector <std::vector <Point <dim> > > actual_points;
1111 typename std::vector <internal::PointValueHistoryImplementation::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
1113 for (; point != point_geometry_data.end (); ++point)
1115 actual_points.push_back (point->support_point_locations);
1117 locations = actual_points;
1125 get_support_locations (locations);
1134 AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
1136 locations = std::vector<Point <dim> > ();
1139 unsigned int n_quadrature_points = quadrature.
size();
1140 std::vector<Point<dim> > evaluation_points;
1143 typename std::vector <internal::PointValueHistoryImplementation::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
1144 for (
unsigned int data_store_index = 0; point != point_geometry_data.end (); ++point, ++data_store_index)
1148 Point <dim> requested_location = point->requested_location;
1150 fe_values.reinit (cell);
1152 evaluation_points = fe_values.get_quadrature_points();
1153 double distance = cell->diameter ();
1154 unsigned int selected_point = 0;
1156 for (
unsigned int q_point = 0; q_point < n_quadrature_points; q_point++)
1158 if (requested_location.
distance (evaluation_points[q_point]) < distance)
1160 selected_point = q_point;
1161 distance = requested_location.
distance (evaluation_points[q_point]);
1165 locations.push_back (evaluation_points[selected_point]);
1174 out <<
"***PointValueHistory status output***\n\n";
1175 out <<
"Closed: " << closed <<
"\n";
1176 out <<
"Cleared: " << cleared <<
"\n";
1177 out <<
"Triangulation_changed: " << triangulation_changed <<
"\n";
1178 out <<
"Have_dof_handler: " << have_dof_handler <<
"\n";
1179 out <<
"Geometric Data" <<
"\n";
1181 typename std::vector <internal::PointValueHistoryImplementation::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
1182 if (point == point_geometry_data.end ())
1184 out <<
"No points stored currently\n";
1190 for (; point != point_geometry_data.end (); ++point)
1192 out <<
"# Requested location: " << point->requested_location <<
"\n";
1193 out <<
"# DoF_index : Support location (for each component)\n";
1194 for (
unsigned int component = 0; component < dof_handler->get_fe(0).n_components (); component++)
1196 out << point->solution_indices[component] <<
" : " << point->support_point_locations [component] <<
"\n";
1203 out <<
"#Cannot access DoF_indices once cleared\n";
1208 if (independent_values.size () != 0)
1210 out <<
"Independent value(s): " << independent_values.size () <<
" : " << independent_values[0].size () <<
"\n";
1211 if (indep_names.size() > 0)
1214 for (
unsigned int name = 0; name < indep_names.size(); name++)
1216 out <<
"<" << indep_names [name] <<
"> ";
1223 out <<
"No independent values stored\n";
1226 std::map <std::string, std::vector <std::vector <double> > >::iterator
1227 data_store_begin = data_store.begin ();
1228 if (data_store_begin != data_store.end())
1230 out <<
"Mnemonic: data set size (mask size, n true components) : n data sets\n";
1232 for (; data_store_begin != data_store.end (); ++data_store_begin)
1235 std::string vector_name = data_store_begin->first;
1236 typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(vector_name);
1237 Assert (mask != component_mask.end(),
ExcMessage(
"vector_name not in class"));
1238 typename std::map <std::string, std::vector <std::string> >::iterator component_names = component_names_map.find(vector_name);
1239 Assert (component_names != component_names_map.end(),
ExcMessage(
"vector_name not in class"));
1241 if (data_store_begin->second.size () != 0)
1243 out << data_store_begin->first <<
": " << data_store_begin->second.size () <<
" (";
1244 out << mask->second.size() <<
", " << mask->second.n_selected_components() <<
") : ";
1245 out << (data_store_begin->second)[0].size () <<
"\n";
1249 out << data_store_begin->first <<
": " << data_store_begin->second.size () <<
" (";
1250 out << mask->second.size() <<
", " << mask->second.n_selected_components() <<
") : ";
1251 out <<
"No points added" <<
"\n";
1254 if (component_names->second.size() > 0)
1256 for (
unsigned int name = 0; name < component_names->second.size(); name++)
1258 out <<
"<" << component_names->second[name] <<
"> ";
1264 out <<
"***end of status output***\n\n";
1279 if (dataset_key.size () != independent_values[0].size ())
1284 std::map <std::string, std::vector <std::vector <double> > >::iterator
1285 data_store_begin = data_store.begin ();
1286 if (have_dof_handler)
1288 for (; data_store_begin != data_store.end (); ++data_store_begin)
1290 Assert (data_store_begin->second.size() > 0,
1292 if ((data_store_begin->second)[0].size () != dataset_key.size ())
1306 if (std::abs ((
int) dataset_key.size () - (int) independent_values[0].size ()) >= 2)
1312 if (have_dof_handler)
1314 std::map <std::string, std::vector <std::vector <double> > >::iterator
1315 data_store_begin = data_store.begin ();
1316 for (; data_store_begin != data_store.end (); ++data_store_begin)
1318 Assert (data_store_begin->second.size() > 0,
1321 if (std::abs ((
int) (data_store_begin->second)[0].size () - (int) dataset_key.size ()) >= 2)
1350 triangulation_changed =
true;
1355 #include "point_value_history.inst" 1358 DEAL_II_NAMESPACE_CLOSE
std::map< std::string, std::vector< std::vector< double > > > data_store
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
void get_postprocessor_locations(const Quadrature< dim > &quadrature, std::vector< Point< dim > > &locations)
PointValueHistory(const unsigned int n_independent_variables=0)
std::vector< std::vector< Tensor< 2, spacedim > > > solution_hessians
virtual void evaluate_scalar_field(const DataPostprocessorInputs::Scalar< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
std::vector< double > solution_values
cell_iterator end() const
const FiniteElement< dim, spacedim > & get_fe() const
std::vector< double > dataset_key
bool represents_the_all_selected_mask() const
Transformed quadrature points.
#define AssertThrow(cond, exc)
std::vector< std::string > indep_names
void add_points(const std::vector< Point< dim > > &locations)
active_cell_iterator begin_active(const unsigned int level=0) const
void add_component_names(const std::string &vector_name, const std::vector< std::string > &component_names)
void add_field_name(const std::string &vector_name, const ComponentMask &component_mask=ComponentMask())
boost::signals2::connection tria_listener
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
std::vector< std::vector< Tensor< 1, spacedim > > > solution_gradients
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void evaluate_field_at_requested_location(const std::string &name, const VectorType &solution)
std::vector< std::vector< double > > independent_values
void evaluate_field(const std::string &name, const VectorType &solution)
void start_new_dataset(const double key)
SmartPointer< const DoFHandler< dim >, PointValueHistory< dim > > dof_handler
PointValueHistory & operator=(const PointValueHistory &point_value_history)
std::vector< Tensor< 2, spacedim > > solution_hessians
Second derivatives of shape functions.
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
unsigned int size() const
void write_gnuplot(const std::string &base_name, const std::vector< Point< dim > > &postprocessor_locations=std::vector< Point< dim > >())
std::vector< Tensor< 1, spacedim > > solution_gradients
void add_independent_names(const std::vector< std::string > &independent_names)
std::vector< internal::PointValueHistoryImplementation::PointGeometryData< dim > > point_geometry_data
void tria_change_listener()
Vector< double > mark_support_locations()
bool deep_check(const bool strict)
std::vector< Point< spacedim > > evaluation_points
bool triangulation_changed
std::map< std::string, ComponentMask > component_mask
Shape function gradients.
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
void status(std::ostream &out)
static ::ExceptionBase & ExcNotImplemented()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void add_point(const Point< dim > &location)
void get_points(std::vector< std::vector< Point< dim > > > &locations)
std::map< std::string, std::vector< std::string > > component_names_map
void get_support_locations(std::vector< std::vector< Point< dim > > > &locations)
virtual std::vector< std::string > get_names() const =0
std::vector<::Vector< double > > solution_values
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)
static ::ExceptionBase & ExcInternalError()
virtual UpdateFlags get_needed_update_flags() const =0
void push_back_independent(const std::vector< double > &independent_values)