adamantine
experimental_data_utils.cc
Go to the documentation of this file.
1 /* SPDX-FileCopyrightText: Copyright (c) 2021 - 2024, the adamantine authors.
2  * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3  */
4 
6 #include <utils.hh>
7 
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>
13 
14 #include <boost/algorithm/string.hpp>
15 
16 #include <fstream>
17 #include <unordered_set>
18 
19 namespace adamantine
20 {
21 template <int dim>
22 std::pair<std::vector<dealii::types::global_dof_index>,
23  std::vector<dealii::Point<dim>>>
24 get_dof_to_support_mapping(dealii::DoFHandler<dim> const &dof_handler)
25 {
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;
29 
30  // Manually do what dealii::DoFTools::map_dofs_to_support_points does, since
31  // that doesn't currently work with FE_Nothing. The issue has been fixed in
32  // deal.II 9.7
33  const dealii::FiniteElement<dim> &fe = dof_handler.get_fe(0);
34 
35  dealii::FEValues<dim, dim> fe_values(fe, fe.get_unit_support_points(),
36  dealii::update_quadrature_points);
37 
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();
41 
42  for (auto const &cell :
43  dealii::filter_iterators(dof_handler.active_cell_iterators(),
44  dealii::IteratorFilters::ActiveFEIndexEqualTo(
45  0, /* locally owned */ true)))
46  {
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)
52  {
53  // Skip duplicate points like vertices and indices that correspond to
54  // ghosted elements
55  if ((visited_dof_indices.count(local_dof_indices[i]) == 0) &&
56  (locally_owned_dofs.is_element(local_dof_indices[i])))
57  {
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]);
61  }
62  }
63  }
64 
65  return {dof_indices, support_points};
66 }
67 
68 template <int dim>
69 std::pair<std::vector<int>, std::vector<int>>
71  dealii::DoFHandler<dim> const &dof_handler)
72 {
73  auto [dof_indices, support_points] = get_dof_to_support_mapping(dof_handler);
74 
75  auto communicator = dof_handler.get_communicator();
76  int my_rank = dealii::Utilities::MPI::this_mpi_process(communicator);
77 
78  // Perform the search
79  dealii::ArborXWrappers::DistributedTree distributed_tree(communicator,
80  support_points);
81  dealii::ArborXWrappers::PointNearestPredicate pt_nearest(points_values.points,
82  1);
83  auto [indices_ranks, offset] = distributed_tree.query(pt_nearest);
84 
85  // Convert the indices and offsets to a pair that maps experimental indices to
86  // dof indices
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)
93  {
94  for (int j = offset[i]; j < offset[i + 1]; ++j)
95  {
96  obs_indices[j] = i;
97  global_indices[j] = indices_ranks[j].second == my_rank
98  ? dof_indices[indices_ranks[j].first]
99  : -1;
100  }
101  }
102 
103  dealii::Utilities::MPI::max(global_indices, communicator, global_indices);
104 
105  return std::make_pair(obs_indices, global_indices);
106 }
107 
108 template <int dim>
110  MPI_Comm const &communicator, PointsValues<dim> const &points_values,
111  std::pair<std::vector<int>, std::vector<int>> &expt_to_dof_mapping,
112  dealii::LinearAlgebra::distributed::Vector<double> &temperature,
113  bool verbose_output)
114 {
115  unsigned int rank = dealii::Utilities::MPI::this_mpi_process(communicator);
116  if ((rank == 0) && verbose_output)
117  {
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())
121  .size()
122  << std::endl;
123  }
124 
125  auto locally_owned_dofs = temperature.locally_owned_elements();
126  for (unsigned int i = 0; i < points_values.values.size(); ++i)
127  {
128  // Only add locally owned elements
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];
131  }
132 
133  temperature.compress(dealii::VectorOperation::insert);
134 }
135 
136 std::vector<std::vector<double>>
137 read_frame_timestamps(boost::property_tree::ptree const &experiment_database)
138 {
139  // PropertyTreeInput experiment.log_filename
140  std::string log_filename =
141  experiment_database.get<std::string>("log_filename");
142 
143  wait_for_file(log_filename, "Waiting for frame time stamps: " + log_filename);
144 
145  // PropertyTreeInput experiment.first_frame_temporal_offset
146  double first_frame_offset =
147  experiment_database.get("first_frame_temporal_offset", 0.0);
148 
149  // PropertyTreeInput experiment.first_frame
150  unsigned int first_frame =
151  experiment_database.get<unsigned int>("first_frame", 0);
152  // PropertyTreeInput experiment.last_frame
153  unsigned int last_frame = experiment_database.get<unsigned int>("last_frame");
154 
155  // PropertyTreeInput experiment.first_camera_id
156  unsigned int first_camera_id =
157  experiment_database.get<unsigned int>("first_camera_id");
158  // PropertyTreeInput experiment.last_camera_id
159  unsigned int last_camera_id =
160  experiment_database.get<unsigned int>("last_camera_id");
161 
162  unsigned int num_cameras = last_camera_id - first_camera_id + 1;
163  std::vector<std::vector<double>> time_stamps(num_cameras);
164 
165  // Read and parse the file
166  std::ifstream file;
167  file.open(log_filename);
168  std::string line;
169  while (std::getline(file, line))
170  {
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())
176  {
177  std::string substring;
178  std::getline(s_stream, substring, ',');
179  boost::trim(substring);
180 
181  if (entry_index == 0)
182  {
183  std::string error_message = "The file " + log_filename +
184  " does not have consecutive frame indices.";
185  ASSERT_THROW(std::stoi(substring) - frame == 1 ||
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;
191  }
192  else
193  {
194  if (frame_of_interest && substring.size() > 0)
195  time_stamps[entry_index - 1].push_back(std::stod(substring) +
196  first_frame_offset);
197  }
198  entry_index++;
199  }
200  }
201  return time_stamps;
202 }
203 
204 } // namespace adamantine
205 
206 //-------------------- Explicit Instantiations --------------------//
207 namespace adamantine
208 {
210  MPI_Comm const &communicator, PointsValues<2> const &points_values,
211  std::pair<std::vector<int>, std::vector<int>> &expt_to_dof_mapping,
212  dealii::LinearAlgebra::distributed::Vector<double> &temperature,
213  bool verbose_output);
215  MPI_Comm const &communicator, PointsValues<3> const &points_values,
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>>>
227 get_dof_to_support_mapping(dealii::DoFHandler<2> const &dof_handler);
228 template std::pair<std::vector<dealii::types::global_dof_index>,
229  std::vector<dealii::Point<3>>>
230 get_dof_to_support_mapping(dealii::DoFHandler<3> const &dof_handler);
231 } // namespace adamantine
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)
Definition: utils.hh:24
void ASSERT_THROW(bool cond, std::string const &message)
Definition: utils.hh:70
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