Flow123d  release_2.2.0-21-g2806cf6
observe.cc
Go to the documentation of this file.
1 /*
2  * observe.cc
3  *
4  * Created on: Jun 28, 2016
5  * Author: jb
6  */
7 
8 #include <string>
9 #include <numeric>
10 #include <cmath>
11 #include <algorithm>
12 #include <unordered_set>
13 
14 #include "system/global_defs.h"
15 #include "input/accessors.hh"
16 #include "input/input_type.hh"
18 
19 #include "mesh/mesh.h"
20 #include "mesh/bih_tree.hh"
21 #include "mesh/region.hh"
22 #include "io/observe.hh"
23 #include "io/element_data_cache.hh"
24 
25 
26 namespace IT = Input::Type;
27 
28 
30  return IT::Record("ObservePoint", "Specification of the observation point. The actual observe element and the observe point on it is determined as follows:\n\n"
31  "1. Find an initial element containing the initial point. If no such element exists we report the error.\n"
32  "2. Use BFS starting from the inital element to find the 'observe element'. The observe element is the closest element "
33  "3. Find the closest projection of the inital point on the observe element and snap this projection according to the 'snap_dim'.\n")
34  .allow_auto_conversion("point")
35  .declare_key("name", IT::String(),
37  "Default name have the form 'obs_<id>', where 'id' "
38  "is the rank of the point on the input."),
39  "Optional point name. Has to be unique. Any string that is valid YAML key in record without any quoting can be used however"
40  "using just alpha-numerical characters and underscore instead of the space is recommended. "
41  )
43  "Initial point for the observe point search.")
44  .declare_key("snap_dim", IT::Integer(0, 4), IT::Default("4"),
45  "The dimension of the sub-element to which center we snap. For value 4 no snapping is done. "
46  "For values 0 up to 3 the element containing the initial point is found and then the observe"
47  "point is snapped to the nearest center of the sub-element of the given dimension. "
48  "E.g. for dimension 2 we snap to the nearest center of the face of the initial element."
49  )
50  .declare_key("snap_region", IT::String(), IT::Default("\"ALL\""),
51  "The region of the initial element for snapping. Without snapping we make a projection to the initial element.")
52  .declare_key("n_search_levels", IT::Integer(0), IT::Default("1"),
53  "Maximum number of levels of the breadth first search used to find the observe element from the initial element. Value zero means to search only the initial element itself.")
54  .close();
55 }
56 
58 {}
59 
60 
61 ObservePoint::ObservePoint(Input::Record in_rec, unsigned int point_idx)
62 : distance_(numeric_limits<double>::infinity())
63 {
64  in_rec_ = in_rec;
65 
66  string default_label = string("obs_") + std::to_string(point_idx);
67  name_ = in_rec.val<string>("name", default_label );
68 
69  vector<double> tmp_coords;
70  in_rec.val<Input::Array>("point").copy_to(tmp_coords);
71  input_point_= arma::vec(tmp_coords);
72 
73  snap_dim_ = in_rec.val<unsigned int>("snap_dim");
74 
75  snap_region_name_ = in_rec.val<string>("snap_region");
76 
77  max_levels_ = in_rec_.val<unsigned int>("n_search_levels");
78 }
79 
80 
81 
82 void ObservePoint::update_projection(unsigned int i_elm, arma::vec local_coords, arma::vec3 global_coords)
83 {
84  double dist = arma::norm(global_coords - input_point_, 2);
85  //cout << "dist: " << dist << endl;
86  if (dist < distance_) {
87  distance_ = dist;
88  element_idx_ = i_elm;
89  local_coords_ = local_coords;
91  }
92 }
93 
94 
95 
97  return distance_ < numeric_limits<double>::infinity();
98 }
99 
100 
101 
102 template <int ele_dim>
104 {
105  if (this->snap_dim_ > ele_dim) return;
106 
107  double min_dist = 2.0; // on the ref element the max distance should be about 1.0, smaler then 2.0
108  arma::vec min_center;
109  for(auto &center : RefElement<ele_dim>::centers_of_subelements(this->snap_dim_))
110  {
111  double dist = arma::norm(center-local_coords_, 2);
112  if ( dist < min_dist) {
113  min_dist = dist;
114  min_center = center;
115  }
116  }
117  this->local_coords_ = min_center;
118 }
119 
120 
122 {
123  Element & elm = mesh.element[element_idx_];
124  switch (elm.dim()) {
125  case 1: snap_to_subelement<1>(); break;
126  case 2: snap_to_subelement<2>(); break;
127  case 3: snap_to_subelement<3>(); break;
128  default: ASSERT(false).error("Clipping supported only for dim=1,2,3.");
129  }
130  this->global_coords_ = elm.element_map() * arma::join_cols(this->local_coords_, arma::ones(1));
131 }
132 
133 
134 
137  if (region_set.size() == 0)
138  THROW( RegionDB::ExcUnknownSet() << RegionDB::EI_Label(snap_region_name_) << in_rec_.ei_address() );
139 
140 
141  const BIHTree &bih_tree=mesh.get_bih_tree();
142  vector<unsigned int> candidate_list, process_list;
143  std::unordered_set<unsigned int> closed_elements(1023);
144 
145  // search for the initial element
146  auto projected_point = bih_tree.tree_box().project_point(input_point_);
147  bih_tree.find_point( projected_point, candidate_list );
148  process_list.swap(candidate_list);
149  candidate_list.clear();
150 
151  unsigned int min_dist_idx=0;
152  double min_dist=numeric_limits<double>::max();
153  for (unsigned int i_candidate=0; i_candidate<process_list.size(); ++i_candidate) {
154  unsigned int i_elm=process_list[i_candidate];
155  Element & elm = mesh.element[i_elm];
156  arma::mat map = elm.element_map();
157 
158  // get barycentric coordinates (1,2,0)
159  arma::vec projection = elm.project_point(input_point_, map);
160 
161  // check that point is on the element
162  if (projection.min() >= -BoundingBox::epsilon) {
163  // This is initial element.
164  //input_point_.print(cout, "input_point");
165  //cout << "i_el: " << i_elm << endl;
166  //projection.print(cout, "projection");
167 
168 
169  // if element match region filter store it as observe element to the obs. point
170  if (elm.region().is_in_region_set(region_set)) {
171  projection[elm.dim()] = 1.0; // use last coordinates for translation
172  arma::vec global_coord = map*projection;
173  update_projection(i_elm, projection.rows(0, elm.dim()-1), global_coord);
174  }
175 
176  closed_elements.insert(i_elm);
177  // add all node neighbours to the next level list
178  for (unsigned int n=0; n < elm.n_nodes(); n++) {
179  for(unsigned int i_node_ele : mesh.node_elements()[mesh.node_vector.index(elm.node[n])])
180  candidate_list.push_back(i_node_ele);
181  }
182  } else {
183  // Point out of the element. Keep the closest element.
184  double distance = fabs(projection.min());
185  if (distance < min_dist) {
186  min_dist=distance;
187  min_dist_idx = i_candidate;
188  }
189  //DebugOut() << print_var(i_candidate);
190  //DebugOut() << print_var(projection);
191  }
192  }
193 
194  if (candidate_list.size() == 0) {
195 
196  unsigned int i_elm=process_list[min_dist_idx];
197  Element & elm = mesh.element[i_elm];
198  // if element match region filter store it as observe element to the obs. point
199  if (elm.region().is_in_region_set(region_set)) {
200  arma::mat map = elm.element_map();
201  arma::vec projection = elm.project_point(input_point_, map);
202  projection = elm.clip_to_element(projection);
203 
204  projection[elm.dim()] = 1.0; // use last coordinates for translation
205  arma::vec global_coord = map*projection;
206  update_projection(i_elm, projection.rows(0, elm.dim()-1), global_coord);
207  }
208 
209  WarningOut().fmt("Failed to find the element containing the initial observe point ({}).\n"
210  "Using the closest element instead.\n", in_rec_.address_string());
211 
212  closed_elements.insert(i_elm);
213  // add all node neighbours to the next level list
214  for (unsigned int n=0; n < elm.n_nodes(); n++) {
215  for(unsigned int i_node_ele : mesh.node_elements()[mesh.node_vector.index(elm.node[n])])
216  candidate_list.push_back(i_node_ele);
217  }
218 
219  }
220 
221  // Try to snap to the observe element with required snap_region
222  for(unsigned int i_level=0; i_level < max_levels_; i_level++) {
223  if (have_observe_element()) break;
224  process_list.swap(candidate_list);
225  candidate_list.clear();
226  for(unsigned int i_elm : process_list) {
227  if (closed_elements.find(i_elm) != closed_elements.end()) continue;
228  Element & elm = mesh.element[i_elm];
229 
230  // if element match region filter, update the obs. point
231  if (elm.region().is_in_region_set(region_set)) {
232  arma::mat map = elm.element_map();
233  arma::vec projection = elm.project_point(input_point_, map);
234  arma::vec point_on_element = elm.clip_to_element(projection);
235 
236  point_on_element[elm.dim()] = 1.0; // use last coordinates for translation
237  arma::vec global_coord = map*point_on_element;
238  update_projection(i_elm, point_on_element.rows(0, elm.dim()-1), global_coord);
239  }
240  // add all node neighbours to the next level list
241  for (unsigned int n=0; n < elm.n_nodes(); n++) {
242  for(unsigned int i_node_ele : mesh.node_elements()[mesh.node_vector.index(elm.node[n])])
243  candidate_list.push_back(i_node_ele);
244  }
245  }
246  }
247  if (! have_observe_element()) {
248  THROW(ExcNoObserveElement() << EI_RegionName(snap_region_name_) << EI_NLevels(max_levels_) );
249  }
250  snap( mesh );
251  Element & elm = mesh.element[element_idx_];
252  double dist = arma::norm(elm.centre() - input_point_, 2);
253  double elm_norm = arma::norm(elm.bounding_box().max() - elm.bounding_box().min(), 2);
254  if (dist > 2*elm_norm)
255  WarningOut().fmt("Observe point ({}) is too distant from the mesh.\n", name_);
256 }
257 
258 
259 
260 void ObservePoint::output(ostream &out, unsigned int indent_spaces, unsigned int precision)
261 {
262  out << setw(indent_spaces) << "" << "- name: " << name_ << endl;
263  out << setw(indent_spaces) << "" << " init_point: " << field_value_to_yaml(input_point_, precision) << endl;
264  out << setw(indent_spaces) << "" << " snap_dim: " << snap_dim_ << endl;
265  out << setw(indent_spaces) << "" << " snap_region: " << snap_region_name_ << endl;
266  out << setw(indent_spaces) << "" << " observe_point: " << field_value_to_yaml(global_coords_, precision) << endl;
267 }
268 
269 
270 
271 
272 Observe::Observe(string observe_name, Mesh &mesh, Input::Array in_array, unsigned int precision)
273 : observe_values_time_(numeric_limits<double>::signaling_NaN()),
274  observe_name_(observe_name),
275  precision_(precision)
276 {
277  // in_rec is Output input record.
278 
279  for(auto it = in_array.begin<Input::Record>(); it != in_array.end(); ++it) {
280  ObservePoint point(*it, points_.size());
281  point.find_observe_point(mesh);
282  points_.push_back( point );
283  observed_element_indices_.push_back(point.element_idx_);
284  }
285  // make indices unique
286  std::sort(observed_element_indices_.begin(), observed_element_indices_.end());
287  auto last = std::unique(observed_element_indices_.begin(), observed_element_indices_.end());
289 
290  time_unit_str_ = "s";
291  time_unit_seconds_ = 1.0;
292 
293  if (points_.size() == 0) return;
295  if (rank_==0) {
296  FilePath observe_file_path(observe_name_ + "_observe.yaml", FilePath::output_file);
297  try {
298  observe_file_path.open_stream(observe_file_);
299  //observe_file_.setf(std::ios::scientific);
300  observe_file_.precision(this->precision_);
301 
302  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, in_array)
303  output_header();
304  }
305 }
306 
308  observe_file_.close();
309 }
310 
311 
312 template <typename T>
313 ElementDataCache<T> & Observe::prepare_compute_data(std::string field_name, double field_time, unsigned int n_rows,
314  unsigned int n_cols)
315 {
317  observe_values_time_ = field_time;
318  else
320  (field_time)(observe_values_time_);
321 
322  OutputDataFieldMap::iterator it=observe_field_values_.find(field_name);
323  if (it == observe_field_values_.end()) {
324  observe_field_values_[field_name]
325  = std::make_shared< ElementDataCache<T> >(field_name, n_rows, n_cols, points_.size());
326  it=observe_field_values_.find(field_name);
327  }
328  return dynamic_cast<ElementDataCache<T> &>(*(it->second));
329 }
330 
331 // explicit instantiation of template method
332 #define OBSERVE_PREPARE_COMPUTE_DATA(TYPE) \
333 template ElementDataCache<TYPE> & Observe::prepare_compute_data<TYPE>(std::string field_name, double field_time, \
334  unsigned int n_rows, unsigned int n_cols)
335 
337 OBSERVE_PREPARE_COMPUTE_DATA(unsigned int);
339 
340 
342  unsigned int indent = 2;
343  observe_file_ << "# Observation file: " << observe_name_ << endl;
344  observe_file_ << "time_unit: " << time_unit_str_ << endl;
345  observe_file_ << "time_unit_in_seconds: " << time_unit_seconds_ << endl;
346  observe_file_ << "points:" << endl;
347  for(auto &point : points_)
348  point.output(observe_file_, indent, precision_);
349  observe_file_ << "data:" << endl;
350 
351 }
352 
353 void Observe::output_time_frame(double time) {
354  if (points_.size() == 0) return;
355 
356  if ( ! no_fields_warning ) {
357  no_fields_warning=true;
358  // check that observe fields are set
360  // first call and no fields
361  ASSERT(observe_field_values_.size() == 0);
362  WarningOut() << "No observe fields for the observation stream: " << observe_name_ << endl;
363  }
364  }
365 
367  ASSERT(observe_field_values_.size() == 0);
368  return;
369  }
370 
371  if (rank_ == 0) {
372  unsigned int indent = 2;
373  observe_file_ << setw(indent) << "" << "- time: " << observe_values_time_ << endl;
374  for(auto &field_data : observe_field_values_) {
375  observe_file_ << setw(indent) << "" << " " << field_data.second->field_input_name() << ": ";
376  field_data.second->print_all_yaml(observe_file_, precision_);
377  observe_file_ << endl;
378  }
379  }
380 
381  observe_values_time_ = numeric_limits<double>::signaling_NaN();
382 
383 }
std::vector< ObservePoint > points_
Full information about observe points.
Definition: observe.hh:219
unsigned int max_levels_
Definition: observe.hh:128
arma::vec project_point(const arma::vec3 &point, const arma::mat &map) const
Definition: elements.cc:205
Iterator< ValueType > begin() const
#define OBSERVE_PREPARE_COMPUTE_DATA(TYPE)
Definition: observe.cc:332
OutputDataFieldMap observe_field_values_
Stored field values.
Definition: observe.hh:224
unsigned int n_nodes() const
vector< vector< unsigned int > > const & node_elements()
Definition: mesh.cc:787
std::ofstream observe_file_
Output file stream.
Definition: observe.hh:234
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
string field_value_to_yaml(const T &mat, unsigned int prec)
arma::vec3 input_point_
Input coordinates of the initial position of the observation point.
Definition: observe.hh:112
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
arma::mat element_map() const
Definition: elements.h:112
double distance_
Definition: observe.hh:141
static const Input::Type::Record & get_input_type()
Definition: observe.cc:29
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:64
bool have_observe_element()
Definition: observe.cc:96
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
Definition: mesh.h:97
void output(ostream &out, unsigned int indent_spaces, unsigned int precision)
Definition: observe.cc:260
Observe(string observe_name, Mesh &mesh, Input::Array in_array, unsigned int precision)
Definition: observe.cc:272
static const double epsilon
stabilization parameter
Definition: bounding_box.hh:57
Node ** node
Definition: elements.h:79
std::string time_unit_str_
String representation of the time unit.
Definition: observe.hh:237
int rank_
Definition: observe.hh:216
const RegionDB & region_db() const
Definition: mesh.h:155
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
std::string to_string(const T &value)
Definition: string.h:29
std::string observe_name_
Definition: observe.hh:231
Class for declaration of the integral input data.
Definition: type_base.hh:489
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:303
Class for declaration of inputs sequences.
Definition: type_base.hh:345
IteratorBase end() const
arma::vec clip_to_element(arma::vec &barycentric)
Definition: elements.h:144
unsigned int dim() const
Input::Record in_rec_
Index in the input array.
Definition: observe.hh:106
EI_Address ei_address() const
Definition: accessors.cc:178
void open_stream(Stream &stream) const
Definition: file_path.cc:211
Global macros to enhance readability and debugging, general constants.
virtual Record & allow_auto_conversion(const string &from_key)
Allows shorter input of the Record providing only value of the from_key given as the parameter...
Definition: type_record.cc:132
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:540
unsigned int precision_
Precision of float output.
Definition: observe.hh:241
void find_observe_point(Mesh &mesh)
Definition: observe.cc:135
ObservePoint()
Definition: observe.cc:57
DummyInt isnan(...)
Definition: format.h:306
void snap_to_subelement()
Definition: observe.cc:103
std::string name_
Observation point name.
Definition: observe.hh:109
Region region() const
Definition: elements.cc:155
const BoundingBox & tree_box() const
Definition: bih_tree.cc:188
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
const Ret val(const string &key) const
void output_header()
Definition: observe.cc:341
std::vector< unsigned int > observed_element_indices_
Elements of the o_points.
Definition: observe.hh:221
const Point & max() const
Class for O(log N) lookup for intersections with a set of bounding boxes.
Definition: bih_tree.hh:36
unsigned int index(const T *pointer) const
Definition: sys_vector.hh:373
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:490
ElementDataCache< T > & prepare_compute_data(std::string field_name, double field_time, unsigned int n_rows, unsigned int n_cols)
Definition: observe.cc:313
double observe_values_time_
Common evaluation time of the fields for single time frame.
Definition: observe.hh:228
void find_point(const Space< 3 >::Point &point, std::vector< unsigned int > &result_list) const
Definition: bih_tree.cc:240
unsigned int snap_dim_
Definition: observe.hh:118
const BIHTree & get_bih_tree()
Definition: mesh.cc:735
#define MPI_Comm_rank
Definition: mpi.h:236
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
static Default read_time(const std::string &description)
The factory function to make an default value that will be specified at the time when a key will be r...
Definition: type_record.hh:97
const double epsilon
Definition: mathfce.h:23
arma::vec3 global_coords_
Global coordinates of the observation point.
Definition: observe.hh:134
arma::vec local_coords_
Local (barycentric) coordinates on the element.
Definition: observe.hh:137
void update_projection(unsigned int i_elm, arma::vec local_coords, arma::vec3 global_coords)
Definition: observe.cc:82
arma::vec3 centre() const
Definition: elements.cc:120
double time_unit_seconds_
Time unit in seconds.
Definition: observe.hh:239
bool is_in_region_set(const RegionSet &set) const
Definition: region.cc:46
Point project_point(const Point &point) const
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:236
#define MPI_COMM_WORLD
Definition: mpi.h:123
const Point & min() const
string snap_region_name_
Definition: observe.hh:123
Record type proxy class.
Definition: type_record.hh:182
unsigned int element_idx_
Final element of the observe point. The index in the mesh.
Definition: observe.hh:131
arma::vec3 global_coords() const
Definition: observe.hh:48
void snap(Mesh &mesh)
Definition: observe.cc:121
bool no_fields_warning
Definition: observe.hh:244
Class for declaration of the input data that are in string format.
Definition: type_base.hh:588
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
~Observe()
Destructor, must close the file.
Definition: observe.cc:307
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:226
void output_time_frame(double time)
Definition: observe.cc:353
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:228
string address_string() const
Definition: accessors.cc:184
BoundingBox bounding_box()
Definition: elements.h:102