8 #include <deal.II/arborx/distributed_tree.h>
9 #include <deal.II/dofs/dof_handler.h>
10 #include <deal.II/fe/mapping_q1.h>
11 #include <deal.II/grid/filtered_iterator.h>
12 #include <deal.II/hp/fe_values.h>
14 #include <boost/algorithm/string.hpp>
17 #include <unordered_set>
22 std::pair<std::vector<dealii::types::global_dof_index>,
23 std::vector<dealii::Point<dim>>>
26 std::vector<dealii::types::global_dof_index> dof_indices;
27 std::vector<dealii::Point<dim>> support_points;
28 std::unordered_set<dealii::types::global_dof_index> visited_dof_indices;
33 const dealii::FiniteElement<dim> &fe = dof_handler.get_fe(0);
35 dealii::FEValues<dim, dim> fe_values(fe, fe.get_unit_support_points(),
36 dealii::update_quadrature_points);
38 std::vector<dealii::types::global_dof_index> local_dof_indices(
39 fe.n_dofs_per_cell());
40 auto locally_owned_dofs = dof_handler.locally_owned_dofs();
42 for (
auto const &cell :
43 dealii::filter_iterators(dof_handler.active_cell_iterators(),
44 dealii::IteratorFilters::ActiveFEIndexEqualTo(
47 fe_values.reinit(cell);
48 cell->get_dof_indices(local_dof_indices);
49 const std::vector<dealii::Point<dim>> &points =
50 fe_values.get_quadrature_points();
51 for (
unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
55 if ((visited_dof_indices.count(local_dof_indices[i]) == 0) &&
56 (locally_owned_dofs.is_element(local_dof_indices[i])))
58 dof_indices.push_back(local_dof_indices[i]);
59 support_points.push_back(points[i]);
60 visited_dof_indices.insert(local_dof_indices[i]);
65 return {dof_indices, support_points};
69 std::pair<std::vector<int>, std::vector<int>>
71 dealii::DoFHandler<dim>
const &dof_handler)
75 auto communicator = dof_handler.get_communicator();
76 int my_rank = dealii::Utilities::MPI::this_mpi_process(communicator);
79 dealii::ArborXWrappers::DistributedTree distributed_tree(communicator,
81 dealii::ArborXWrappers::PointNearestPredicate pt_nearest(points_values.
points,
83 auto [indices_ranks, offset] = distributed_tree.query(pt_nearest);
87 std::vector<int> obs_indices;
88 std::vector<int> global_indices;
89 obs_indices.resize(indices_ranks.size());
90 global_indices.resize(indices_ranks.size());
91 unsigned int const n_queries = points_values.
points.size();
92 for (
unsigned int i = 0; i < n_queries; ++i)
94 for (
int j = offset[i]; j < offset[i + 1]; ++j)
97 global_indices[j] = indices_ranks[j].second == my_rank
98 ? dof_indices[indices_ranks[j].first]
103 dealii::Utilities::MPI::max(global_indices, communicator, global_indices);
105 return std::make_pair(obs_indices, global_indices);
111 std::pair<std::vector<int>, std::vector<int>> &expt_to_dof_mapping,
112 dealii::LinearAlgebra::distributed::Vector<double> &temperature,
115 unsigned int rank = dealii::Utilities::MPI::this_mpi_process(communicator);
116 if ((rank == 0) && verbose_output)
118 std::cout <<
"Number of unique experimental points is "
119 << std::set<double>(expt_to_dof_mapping.second.begin(),
120 expt_to_dof_mapping.second.end())
125 auto locally_owned_dofs = temperature.locally_owned_elements();
126 for (
unsigned int i = 0; i < points_values.
values.size(); ++i)
129 if (locally_owned_dofs.is_element(expt_to_dof_mapping.second[i]))
130 temperature[expt_to_dof_mapping.second[i]] = points_values.
values[i];
133 temperature.compress(dealii::VectorOperation::insert);
136 std::vector<std::vector<double>>
140 std::string log_filename =
141 experiment_database.get<std::string>(
"log_filename");
143 wait_for_file(log_filename,
"Waiting for frame time stamps: " + log_filename);
146 double first_frame_offset =
147 experiment_database.get(
"first_frame_temporal_offset", 0.0);
150 unsigned int first_frame =
151 experiment_database.get<
unsigned int>(
"first_frame", 0);
153 unsigned int last_frame = experiment_database.get<
unsigned int>(
"last_frame");
156 unsigned int first_camera_id =
157 experiment_database.get<
unsigned int>(
"first_camera_id");
159 unsigned int last_camera_id =
160 experiment_database.get<
unsigned int>(
"last_camera_id");
162 unsigned int num_cameras = last_camera_id - first_camera_id + 1;
163 std::vector<std::vector<double>> time_stamps(num_cameras);
167 file.open(log_filename);
169 while (std::getline(file,
line))
171 unsigned int entry_index = 0;
172 std::stringstream s_stream(
line);
173 bool frame_of_interest =
false;
174 unsigned int frame = std::numeric_limits<unsigned int>::max();
175 while (s_stream.good())
177 std::string substring;
178 std::getline(s_stream, substring,
',');
179 boost::trim(substring);
181 if (entry_index == 0)
183 std::string error_message =
"The file " + log_filename +
184 " does not have consecutive frame indices.";
186 frame == std::numeric_limits<unsigned int>::max(),
187 error_message.c_str());
188 frame = std::stoi(substring);
189 if (frame >= first_frame && frame <= last_frame)
190 frame_of_interest =
true;
194 if (frame_of_interest && substring.size() > 0)
195 time_stamps[entry_index - 1].push_back(std::stod(substring) +
211 std::pair<std::vector<int>, std::vector<int>> &expt_to_dof_mapping,
212 dealii::LinearAlgebra::distributed::Vector<double> &temperature,
213 bool verbose_output);
216 std::pair<std::vector<int>, std::vector<int>> &expt_to_dof_mapping,
217 dealii::LinearAlgebra::distributed::Vector<double> &temperature,
218 bool verbose_output);
219 template std::pair<std::vector<int>, std::vector<int>>
221 dealii::DoFHandler<2>
const &dof_handler);
222 template std::pair<std::vector<int>, std::vector<int>>
224 dealii::DoFHandler<3>
const &dof_handler);
225 template std::pair<std::vector<dealii::types::global_dof_index>,
226 std::vector<dealii::Point<2>>>
228 template std::pair<std::vector<dealii::types::global_dof_index>,
229 std::vector<dealii::Point<3>>>
std::pair< std::vector< dealii::types::global_dof_index >, std::vector< dealii::Point< dim > > > get_dof_to_support_mapping(dealii::DoFHandler< dim > const &dof_handler)
std::pair< std::vector< int >, std::vector< int > > get_expt_to_dof_mapping(PointsValues< dim > const &points_values, dealii::DoFHandler< dim > const &dof_handler)
void wait_for_file(std::string const &filename, std::string const &message)
void ASSERT_THROW(bool cond, std::string const &message)
void set_with_experimental_data(MPI_Comm const &communicator, PointsValues< dim > const &points_values, std::pair< std::vector< int >, std::vector< int >> &expt_to_dof_mapping, dealii::LinearAlgebra::distributed::Vector< double > &temperature, bool verbose_output)
std::vector< std::vector< double > > read_frame_timestamps(boost::property_tree::ptree const &experiment_database)
std::vector< dealii::Point< dim > > points
std::vector< double > values