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