Flow123d  release_2.2.0-914-gf1a3a4f
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 #include <queue>
14 
15 #include "system/global_defs.h"
16 #include "input/accessors.hh"
17 #include "input/input_type.hh"
19 
20 #include "mesh/mesh.h"
21 #include "mesh/bih_tree.hh"
22 #include "mesh/region.hh"
23 #include "io/observe.hh"
24 #include "io/element_data_cache.hh"
25 #include "fem/mapping_p1.hh"
26 
27 
28 namespace IT = Input::Type;
29 
30 
31 /**
32  * Helper class allows to work with ObservePoint above elements with different dimensions
33  *
34  * Allows:
35  * - calculate projection of points by dimension
36  * - snap to subelement
37  */
38 template<unsigned int dim>
40 public:
41  /// Constructor
43 
44  ObservePointData projection(arma::vec3 input_point, unsigned int i_elm, Element &elm) {
45  arma::mat::fixed<3, dim+1> elm_map = mapping_.element_map(elm);
46  arma::vec::fixed<dim+1> projection = mapping_.project_real_to_unit(input_point, elm_map);
47  projection = mapping_.clip_to_element(projection);
48 
49  ObservePointData data;
50  data.element_idx_ = i_elm;
51  data.local_coords_ = projection.rows(1, elm.dim());
52  data.global_coords_ = elm_map*projection;
53  data.distance_ = arma::norm(data.global_coords_ - input_point, 2);
54 
55  return data;
56  }
57 
58  /**
59  * Snap local coords to the subelement. Called by the ObservePoint::snap method.
60  */
61  void snap_to_subelement(ObservePointData & observe_data, Element & elm, unsigned int snap_dim)
62  {
63  if (snap_dim <= dim) {
64  double min_dist = 2.0; // on the ref element the max distance should be about 1.0, smaler then 2.0
65  arma::vec min_center;
66  for(auto &center : RefElement<dim>::centers_of_subelements(snap_dim))
67  {
68  double dist = arma::norm(center - observe_data.local_coords_, 2);
69  if ( dist < min_dist) {
70  min_dist = dist;
71  min_center = center;
72  }
73  }
74  observe_data.local_coords_ = min_center;
75  }
76 
77  arma::mat::fixed<3, dim+1> elm_map = mapping_.element_map(elm);
78  observe_data.global_coords_ = elm_map * RefElement<dim>::local_to_bary(observe_data.local_coords_);
79  }
80 
81 private:
82  /// Mapping object.
84 };
85 
86 template class ProjectionHandler<1>;
87 template class ProjectionHandler<2>;
88 template class ProjectionHandler<3>;
89 
90 
91 /**
92  * Helper struct, used as comparator of priority queue in ObservePoint::find_observe_point.
93  */
95 {
96  bool operator()(const ObservePointData& lhs, const ObservePointData& rhs) const
97  {
98  return lhs.distance_ > rhs.distance_;
99  }
100 };
101 
102 
103 /*******************************************************************
104  * implementation of ObservePoint
105  */
106 
108  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"
109  "1. Find an initial element containing the initial point. If no such element exists we report the error.\n"
110  "2. Use BFS starting from the inital element to find the 'observe element'. The observe element is the closest element "
111  "3. Find the closest projection of the inital point on the observe element and snap this projection according to the 'snap_dim'.\n")
112  .allow_auto_conversion("point")
113  .declare_key("name", IT::String(),
115  "Default name have the form 'obs_<id>', where 'id' "
116  "is the rank of the point on the input."),
117  "Optional point name. Has to be unique. Any string that is valid YAML key in record without any quoting can be used however"
118  "using just alpha-numerical characters and underscore instead of the space is recommended. "
119  )
121  "Initial point for the observe point search.")
122  .declare_key("snap_dim", IT::Integer(0, 4), IT::Default("4"),
123  "The dimension of the sub-element to which center we snap. For value 4 no snapping is done. "
124  "For values 0 up to 3 the element containing the initial point is found and then the observe"
125  "point is snapped to the nearest center of the sub-element of the given dimension. "
126  "E.g. for dimension 2 we snap to the nearest center of the face of the initial element."
127  )
128  .declare_key("snap_region", IT::String(), IT::Default("\"ALL\""),
129  "The region of the initial element for snapping. Without snapping we make a projection to the initial element.")
130  .declare_key("search_radius", IT::Double(0.0),
131  IT::Default::read_time("Maximal distance of observe point from Mesh relative to its size (bounding box). "),
132  "Global value is define in Mesh by the key global_observe_search_radius.")
133  .close();
134 }
135 
137 {}
138 
139 
140 ObservePoint::ObservePoint(Input::Record in_rec, Mesh &mesh, unsigned int point_idx)
141 {
142  in_rec_ = in_rec;
143 
144  string default_label = string("obs_") + std::to_string(point_idx);
145  name_ = in_rec.val<string>("name", default_label );
146 
147  vector<double> tmp_coords;
148  in_rec.val<Input::Array>("point").copy_to(tmp_coords);
149  input_point_= arma::vec(tmp_coords);
150 
151  snap_dim_ = in_rec.val<unsigned int>("snap_dim");
152 
153  snap_region_name_ = in_rec.val<string>("snap_region");
154 
155  const BoundingBox &main_box = mesh.get_bih_tree().tree_box();
156  double max_mesh_size = arma::max(main_box.max() - main_box.min());
157  max_search_radius_ = in_rec_.val<double>("search_radius", mesh.global_observe_radius()) * max_mesh_size;
158 }
159 
160 
161 
163  return observe_data_.distance_ < numeric_limits<double>::infinity();
164 }
165 
166 
167 
169 {
170  Element & elm = mesh.element[observe_data_.element_idx_];
171  switch (elm.dim()) {
172  case 1:
173  {
175  ph.snap_to_subelement(observe_data_, elm, snap_dim_);
176  break;
177  }
178  case 2:
179  {
181  ph.snap_to_subelement(observe_data_, elm, snap_dim_);
182  break;
183  }
184  case 3:
185  {
187  ph.snap_to_subelement(observe_data_, elm, snap_dim_);
188  break;
189  }
190  default: ASSERT(false).error("Clipping supported only for dim=1,2,3.");
191  }
192 }
193 
194 
195 
197  RegionSet region_set = mesh.region_db().get_region_set(snap_region_name_);
198  if (region_set.size() == 0)
199  THROW( RegionDB::ExcUnknownSet() << RegionDB::EI_Label(snap_region_name_) << in_rec_.ei_address() );
200 
201 
202  const BIHTree &bih_tree=mesh.get_bih_tree();
203  vector<unsigned int> candidate_list;
204  std::unordered_set<unsigned int> closed_elements(1023);
205  std::priority_queue< ObservePointData, std::vector<ObservePointData>, CompareByDist > candidate_queue;
206 
207  // search for the initial element
208  auto projected_point = bih_tree.tree_box().project_point(input_point_);
209  bih_tree.find_point(projected_point, candidate_list, true);
210 
211  for (unsigned int i_candidate=0; i_candidate<candidate_list.size(); ++i_candidate) {
212  unsigned int i_elm=candidate_list[i_candidate];
213  Element & elm = mesh.element[i_elm];
214 
215  // project point, add candidate to queue
216  auto observe_data = point_projection(i_elm, elm);
217  if (observe_data.distance_ <= max_search_radius_)
218  candidate_queue.push(observe_data);
219  closed_elements.insert(i_elm);
220  }
221 
222  while (!candidate_queue.empty())
223  {
224  auto candidate_data = candidate_queue.top();
225  candidate_queue.pop();
226 
227  unsigned int i_elm=candidate_data.element_idx_;
228  Element & elm = mesh.element[i_elm];
229 
230  // test if candidate is in region and update projection
231  if (elm.region().is_in_region_set(region_set)) {
232  ASSERT_LE(candidate_data.distance_, observe_data_.distance_).error();
233 
234  observe_data_.distance_ = candidate_data.distance_;
235  observe_data_.element_idx_ = candidate_data.element_idx_;
236  observe_data_.local_coords_ = candidate_data.local_coords_;
237  observe_data_.global_coords_ = candidate_data.global_coords_;
238  break;
239  }
240 
241  // add candidates to queue
242  for (unsigned int n=0; n < elm.n_nodes(); n++)
243  for(unsigned int i_node_ele : mesh.node_elements()[mesh.node_vector.index(elm.node[n])]) {
244  if (closed_elements.find(i_node_ele) == closed_elements.end()) {
245  Element & neighbor_elm = mesh.element[i_node_ele];
246  auto observe_data = point_projection(i_node_ele, neighbor_elm);
247  if (observe_data.distance_ <= max_search_radius_)
248  candidate_queue.push(observe_data);
249  closed_elements.insert(i_node_ele);
250  }
251  }
252  }
253 
254  if (! have_observe_element()) {
255  THROW(ExcNoObserveElement() << EI_RegionName(snap_region_name_) );
256  }
257  snap( mesh );
258  Element & elm = mesh.element[observe_data_.element_idx_];
259  double dist = arma::norm(elm.centre() - input_point_, 2);
260  double elm_norm = arma::norm(elm.bounding_box().max() - elm.bounding_box().min(), 2);
261  if (dist > 2*elm_norm)
262  WarningOut().fmt("Observe point ({}) is too distant from the mesh.\n", name_);
263 }
264 
265 
266 
267 
268 void ObservePoint::output(ostream &out, unsigned int indent_spaces, unsigned int precision)
269 {
270  out << setw(indent_spaces) << "" << "- name: " << name_ << endl;
271  out << setw(indent_spaces) << "" << " init_point: " << field_value_to_yaml(input_point_) << endl;
272  out << setw(indent_spaces) << "" << " snap_dim: " << snap_dim_ << endl;
273  out << setw(indent_spaces) << "" << " snap_region: " << snap_region_name_ << endl;
274  out << setw(indent_spaces) << "" << " observe_point: " << field_value_to_yaml(observe_data_.global_coords_, precision) << endl;
275 }
276 
277 
278 
280  switch (elm.dim()) {
281  case 1:
282  {
284  return ph.projection(input_point_, i_elm, elm);
285  break;
286  }
287  case 2:
288  {
290  return ph.projection(input_point_, i_elm, elm);
291  break;
292  }
293  case 3:
294  {
296  return ph.projection(input_point_, i_elm, elm);
297  break;
298  }
299  default:
300  ASSERT(false).error("Invalid element dimension!");
301  }
302 
303  return ObservePointData(); // Should not happen.
304 }
305 
306 
307 
308 
309 /*******************************************************************
310  * implementation of Observe
311  */
312 
313 Observe::Observe(string observe_name, Mesh &mesh, Input::Array in_array, unsigned int precision)
314 : observe_values_time_(numeric_limits<double>::signaling_NaN()),
315  observe_name_(observe_name),
316  precision_(precision)
317 {
318  // in_rec is Output input record.
319 
320  for(auto it = in_array.begin<Input::Record>(); it != in_array.end(); ++it) {
321  ObservePoint point(*it, mesh, points_.size());
322  point.find_observe_point(mesh);
323  points_.push_back( point );
324  observed_element_indices_.push_back(point.observe_data_.element_idx_);
325  }
326  // make indices unique
327  std::sort(observed_element_indices_.begin(), observed_element_indices_.end());
328  auto last = std::unique(observed_element_indices_.begin(), observed_element_indices_.end());
330 
331  time_unit_str_ = "s";
332  time_unit_seconds_ = 1.0;
333 
334  if (points_.size() == 0) return;
336  if (rank_==0) {
337  FilePath observe_file_path(observe_name_ + "_observe.yaml", FilePath::output_file);
338  try {
339  observe_file_path.open_stream(observe_file_);
340  //observe_file_.setf(std::ios::scientific);
341  observe_file_.precision(this->precision_);
342 
343  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, in_array)
344  output_header();
345  }
346 }
347 
349  observe_file_.close();
350 }
351 
352 
353 template <typename T>
354 ElementDataCache<T> & Observe::prepare_compute_data(std::string field_name, double field_time, unsigned int n_rows,
355  unsigned int n_cols)
356 {
358  observe_values_time_ = field_time;
359  else
361  (field_time)(observe_values_time_);
362 
363  OutputDataFieldMap::iterator it=observe_field_values_.find(field_name);
364  if (it == observe_field_values_.end()) {
365  observe_field_values_[field_name]
366  = std::make_shared< ElementDataCache<T> >(field_name, n_rows, n_cols, points_.size());
367  it=observe_field_values_.find(field_name);
368  }
369  return dynamic_cast<ElementDataCache<T> &>(*(it->second));
370 }
371 
372 // explicit instantiation of template method
373 #define OBSERVE_PREPARE_COMPUTE_DATA(TYPE) \
374 template ElementDataCache<TYPE> & Observe::prepare_compute_data<TYPE>(std::string field_name, double field_time, \
375  unsigned int n_rows, unsigned int n_cols)
376 
378 OBSERVE_PREPARE_COMPUTE_DATA(unsigned int);
380 
381 
383  unsigned int indent = 2;
384  observe_file_ << "# Observation file: " << observe_name_ << endl;
385  observe_file_ << "time_unit: " << time_unit_str_ << endl;
386  observe_file_ << "time_unit_in_seconds: " << time_unit_seconds_ << endl;
387  observe_file_ << "points:" << endl;
388  for(auto &point : points_)
389  point.output(observe_file_, indent, precision_);
390  observe_file_ << "data:" << endl;
391 
392 }
393 
394 void Observe::output_time_frame(double time) {
395  if (points_.size() == 0) return;
396 
397  if ( ! no_fields_warning ) {
398  no_fields_warning=true;
399  // check that observe fields are set
401  // first call and no fields
402  ASSERT(observe_field_values_.size() == 0);
403  WarningOut() << "No observe fields for the observation stream: " << observe_name_ << endl;
404  }
405  }
406 
408  ASSERT(observe_field_values_.size() == 0);
409  return;
410  }
411 
412  if (rank_ == 0) {
413  unsigned int indent = 2;
414  observe_file_ << setw(indent) << "" << "- time: " << observe_values_time_ << endl;
415  for(auto &field_data : observe_field_values_) {
416  observe_file_ << setw(indent) << "" << " " << field_data.second->field_input_name() << ": ";
417  field_data.second->print_all_yaml(observe_file_, precision_);
418  observe_file_ << endl;
419  }
420  }
421 
422  observe_values_time_ = numeric_limits<double>::signaling_NaN();
423 
424 }
std::vector< ObservePoint > points_
Full information about observe points.
Definition: observe.hh:226
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Iterator< ValueType > begin() const
Bounding box in 3d ambient space.
Definition: bounding_box.hh:45
#define OBSERVE_PREPARE_COMPUTE_DATA(TYPE)
Definition: observe.cc:373
OutputDataFieldMap observe_field_values_
Stored field values.
Definition: observe.hh:231
unsigned int n_nodes() const
vector< vector< unsigned int > > const & node_elements()
Definition: mesh.cc:814
std::ofstream observe_file_
Output file stream.
Definition: observe.hh:241
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
string field_value_to_yaml(const T &mat, unsigned int prec)
BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map) const
Definition: mapping_p1.cc:275
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
ProjectionHandler()
Constructor.
Definition: observe.cc:42
static const Input::Type::Record & get_input_type()
Definition: observe.cc:107
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:64
bool have_observe_element()
Definition: observe.cc:162
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:99
void output(ostream &out, unsigned int indent_spaces, unsigned int precision)
Definition: observe.cc:268
Observe(string observe_name, Mesh &mesh, Input::Array in_array, unsigned int precision)
Definition: observe.cc:313
#define ASSERT_LE(a, b)
Definition of comparative assert macro (Less or Equal)
Definition: asserts.hh:303
void snap_to_subelement(ObservePointData &observe_data, Element &elm, unsigned int snap_dim)
Definition: observe.cc:61
Node ** node
Definition: elements.h:90
unsigned int element_idx_
Final element of the observe point. The index in the mesh.
Definition: observe.hh:30
std::string time_unit_str_
String representation of the time unit.
Definition: observe.hh:244
int rank_
Definition: observe.hh:223
const RegionDB & region_db() const
Definition: mesh.h:170
#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:238
Class for declaration of the integral input data.
Definition: type_base.hh:489
ObservePointData point_projection(unsigned int i_elm, Element &elm)
Project point to given element by dimension of this element.
Definition: observe.cc:279
arma::vec local_coords_
Local (barycentric) coordinates on the element.
Definition: observe.hh:39
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
unsigned int dim() const
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
arma::vec3 global_coords_
Global coordinates of the observation point.
Definition: observe.hh:36
unsigned int precision_
Precision of float output.
Definition: observe.hh:248
void find_observe_point(Mesh &mesh)
Definition: observe.cc:196
DummyInt isnan(...)
Definition: format.h:306
ObservePointData projection(arma::vec3 input_point, unsigned int i_elm, Element &elm)
Definition: observe.cc:44
Region region() const
Definition: elements.cc:165
bool operator()(const ObservePointData &lhs, const ObservePointData &rhs) const
Definition: observe.cc:96
const BoundingBox & tree_box() const
Definition: bih_tree.cc:210
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:382
std::vector< unsigned int > observed_element_indices_
Elements of the o_points.
Definition: observe.hh:228
const Point & max() const
double distance_
Definition: observe.hh:43
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:354
double observe_values_time_
Common evaluation time of the fields for single time frame.
Definition: observe.hh:235
ElementMap element_map(const Element &elm) const
Definition: mapping_p1.cc:265
const BIHTree & get_bih_tree()
Getter for BIH. Creates and compute BIH at first call.
Definition: mesh.cc:752
#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
void find_point(const Space< 3 >::Point &point, std::vector< unsigned int > &result_list, bool full_list=false) const
Definition: bih_tree.cc:265
arma::vec3 centre() const
Computes the barycenter.
Definition: elements.cc:130
double time_unit_seconds_
Time unit in seconds.
Definition: observe.hh:246
bool is_in_region_set(const RegionSet &set) const
Definition: region.cc:46
MappingP1< dim, 3 > mapping_
Mapping object.
Definition: observe.cc:83
Point project_point(const Point &point) const
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
#define MPI_COMM_WORLD
Definition: mpi.h:123
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
Definition: ref_element.cc:286
const Point & min() const
Record type proxy class.
Definition: type_record.hh:182
BaryPoint clip_to_element(BaryPoint &barycentric)
Definition: mapping_p1.cc:297
void snap(Mesh &mesh)
Definition: observe.cc:168
bool no_fields_warning
Definition: observe.hh:251
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:348
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:258
double global_observe_radius() const
Maximal distance of observe point from Mesh relative to its size.
Definition: mesh.cc:758
void output_time_frame(double time)
Definition: observe.cc:394
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:260
BoundingBox bounding_box()
Definition: elements.h:117