Flow123d  last_with_con_2.0.0-4-g42e6930
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/output_data.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;
90  global_coords_ = global_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  bih_tree.find_point(input_point_, candidate_list);
147  process_list.swap(candidate_list);
148  candidate_list.clear();
149  for(unsigned int i_elm : process_list) {
150  Element & elm = mesh.element[i_elm];
151  arma::mat map = elm.element_map();
152  arma::vec projection = elm.project_point(input_point_, map);
153 
154  // check that point is on the element
155  if (projection.min() >=0.0) {
156  // This is initial element.
157  //input_point_.print(cout, "input_point");
158  //cout << "i_el: " << i_elm << endl;
159  //projection.print(cout, "projection");
160 
161  // if element match region filter store it as observe element to the obs. point
162  if (elm.region().is_in_region_set(region_set)) {
163  projection[elm.dim()] = 1.0; // use last coordinates for translation
164  arma::vec global_coord = map*projection;
165  update_projection(i_elm, projection.rows(0, elm.dim()-1), global_coord);
166  }
167 
168  closed_elements.insert(i_elm);
169  // add all node neighbours to the next level list
170  for (unsigned int n=0; n < elm.n_nodes(); n++) {
171  for(unsigned int i_node_ele : mesh.node_elements[mesh.node_vector.index(elm.node[n])])
172  candidate_list.push_back(i_node_ele);
173  }
174  }
175  }
176 
177  if (candidate_list.size() == 0) THROW( ExcNoInitialPoint() << in_rec_.ei_address() );
178 
179  // Try to snap to the observe element with required snap_region
180  for(unsigned int i_level=0; i_level < max_levels_; i_level++) {
181  if (have_observe_element()) break;
182  process_list.swap(candidate_list);
183  candidate_list.clear();
184  for(unsigned int i_elm : process_list) {
185  if (closed_elements.find(i_elm) != closed_elements.end()) continue;
186  Element & elm = mesh.element[i_elm];
187 
188  // if element match region filter, update the obs. point
189  if (elm.region().is_in_region_set(region_set)) {
190  arma::mat map = elm.element_map();
191  arma::vec projection = elm.project_point(input_point_, map);
192  arma::vec point_on_element = elm.clip_to_element(projection);
193 
194  point_on_element[elm.dim()] = 1.0; // use last coordinates for translation
195  arma::vec global_coord = map*point_on_element;
196  update_projection(i_elm, point_on_element.rows(0, elm.dim()-1), global_coord);
197  }
198  // add all node neighbours to the next level list
199  for (unsigned int n=0; n < elm.n_nodes(); n++) {
200  for(unsigned int i_node_ele : mesh.node_elements[mesh.node_vector.index(elm.node[n])])
201  candidate_list.push_back(i_node_ele);
202  }
203  }
204  }
205  if (! have_observe_element()) {
206  THROW(ExcNoObserveElement() << EI_RegionName(snap_region_name_) << EI_NLevels(max_levels_) );
207  }
208  snap( mesh );
209 }
210 
211 
212 
213 void ObservePoint::output(ostream &out, unsigned int indent_spaces, unsigned int precision)
214 {
215  out << setw(indent_spaces) << "" << "- name: " << name_ << endl;
216  out << setw(indent_spaces) << "" << " init_point: " << field_value_to_yaml(input_point_) << endl;
217  out << setw(indent_spaces) << "" << " snap_dim: " << snap_dim_ << endl;
218  out << setw(indent_spaces) << "" << " snap_region: " << snap_region_name_ << endl;
219  out << setw(indent_spaces) << "" << " observe_point: " << field_value_to_yaml(global_coords_, precision) << endl;
220 }
221 
222 
223 
224 
225 Observe::Observe(string observe_name, Mesh &mesh, Input::Array in_array, unsigned int precision)
226 : mesh_(&mesh),
227  observe_values_time_(numeric_limits<double>::signaling_NaN()),
228  precision_(precision)
229 {
230  // in_rec is Output input record.
231 
232  for(auto it = in_array.begin<Input::Record>(); it != in_array.end(); ++it) {
233  ObservePoint point(*it, points_.size());
234  point.find_observe_point(*mesh_);
235  points_.push_back( point );
236  observed_element_indices_.push_back(point.element_idx_);
237  }
238  // make indices unique
239  std::sort(observed_element_indices_.begin(), observed_element_indices_.end());
240  auto last = std::unique(observed_element_indices_.begin(), observed_element_indices_.end());
242 
243  time_unit_str_ = "s";
244  time_unit_seconds_ = 1.0;
245 
246  if (points_.size() == 0) return;
248  if (rank_==0) {
249  FilePath observe_file_path(observe_name + "_observe.yaml", FilePath::output_file);
250  observe_file_path.open_stream(observe_file_);
251  output_header(observe_name);
252  }
253 }
254 
256  observe_file_.close();
257 }
258 
259 
260 template<int spacedim, class Value>
262 {
263  if (points_.size() == 0) return;
264 
265  double field_time = field.time();
267  observe_values_time_ = field_time;
268  else
270  (field_time)(observe_values_time_);
271 
272  OutputDataFieldMap::iterator it=observe_field_values_.find(field.name());
273  if (it == observe_field_values_.end()) {
274  observe_field_values_[field.name()] = std::make_shared< OutputData<Value> >(field, points_.size());
275  it=observe_field_values_.find(field.name());
276  }
277  OutputData<Value> &output_data = dynamic_cast<OutputData<Value> &>(*(it->second));
278 
279  unsigned int i_data=0;
280  for(ObservePoint &o_point : points_) {
281  unsigned int ele_index = o_point.element_idx_;
282  const Value &obs_value =
283  Value( const_cast<typename Value::return_type &>(
284  field.value(o_point.global_coords_,
285  ElementAccessor<spacedim>(this->mesh_, ele_index,false)) ));
286  output_data.store_value(i_data, obs_value);
287  i_data++;
288  }
289 
290 }
291 
292 // Instantiation of the method template for particular dimension.
293 #define INSTANCE_DIM(dim) \
294 template void Observe::compute_field_values(Field<dim, FieldValue<0>::Enum> &); \
295 template void Observe::compute_field_values(Field<dim, FieldValue<0>::Integer> &); \
296 template void Observe::compute_field_values(Field<dim, FieldValue<0>::Scalar> &); \
297 template void Observe::compute_field_values(Field<dim, FieldValue<dim>::VectorFixed> &); \
298 template void Observe::compute_field_values(Field<dim, FieldValue<dim>::TensorFixed> &);
299 
300 // Make all instances for both dimensions.
301 INSTANCE_DIM(2)
302 INSTANCE_DIM(3)
303 
304 
305 void Observe::output_header(string observe_name) {
306  unsigned int indent = 2;
307  observe_file_ << "# Observation file: " << observe_name << endl;
308  observe_file_ << "time_unit: " << time_unit_str_ << endl;
309  observe_file_ << "time_unit_in_secodns: " << time_unit_seconds_ << endl;
310  observe_file_ << "points:" << endl;
311  for(auto &point : points_)
312  point.output(observe_file_, indent, precision_);
313  observe_file_ << "data:" << endl;
314 
315 }
316 
317 void Observe::output_time_frame(double time) {
318  if (points_.size() == 0) return;
319 
320  if (rank_ == 0) {
321  unsigned int indent = 2;
322  observe_file_ << setw(indent) << "" << "- time: " << observe_values_time_ << endl;
323  for(auto &field_data : observe_field_values_) {
324  observe_file_ << setw(indent) << "" << " " << field_data.second->field_name << ": ";
325  field_data.second->print_all_yaml(observe_file_, precision_);
326  observe_file_ << endl;
327  }
328  }
329 
330  observe_values_time_ = numeric_limits<double>::signaling_NaN();
331 
332 }
std::vector< ObservePoint > points_
Full information about observe points.
Definition: observe.hh:192
unsigned int max_levels_
Definition: observe.hh:116
arma::vec project_point(const arma::vec3 &point, const arma::mat &map) const
Definition: elements.cc:205
Iterator< ValueType > begin() const
OutputDataFieldMap observe_field_values_
Stored field values.
Definition: observe.hh:200
unsigned int n_nodes() const
std::ofstream observe_file_
Output file stream.
Definition: observe.hh:207
vector< vector< unsigned int > > node_elements
Definition: mesh.h:280
Accessor to input data conforming to declared Array.
Definition: accessors.hh:552
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:100
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:50
arma::mat element_map() const
Definition: elements.h:112
double distance_
Definition: observe.hh:129
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:62
static const Input::Type::Record & get_input_type()
Definition: observe.cc:29
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:99
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
Definition: mesh.h:95
void output(ostream &out, unsigned int indent_spaces, unsigned int precision)
Definition: observe.cc:213
Observe(string observe_name, Mesh &mesh, Input::Array in_array, unsigned int precision)
Definition: observe.cc:225
Node ** node
Definition: elements.h:79
std::string time_unit_str_
String representation of the time unit.
Definition: observe.hh:210
int rank_
Definition: observe.hh:186
const RegionDB & region_db() const
Definition: mesh.h:152
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
Mesh * mesh_
Definition: observe.hh:189
This class is used for storing data that are copied from field.
Definition: output_data.hh:23
std::string to_string(const T &value)
Definition: string.h:29
Class for declaration of the integral input data.
Definition: type_base.hh:465
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:286
Class for declaration of inputs sequences.
Definition: type_base.hh:321
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:94
EI_Address ei_address() const
Definition: accessors.cc:170
void open_stream(Stream &stream) const
Definition: file_path.cc:211
double time() const
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:518
unsigned int precision_
Precision of float output.
Definition: observe.hh:214
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:97
Region region() const
Definition: elements.cc:155
Accessor to the data with type Type::Record.
Definition: accessors.hh:277
const Ret val(const string &key) const
std::vector< unsigned int > observed_element_indices_
Elements of the o_points.
Definition: observe.hh:194
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
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:342
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:468
double observe_values_time_
Time of fields when the observe values were computed.
Definition: observe.hh:204
void find_point(const Space< 3 >::Point &point, std::vector< unsigned int > &result_list) const
Definition: bih_tree.cc:240
void output_header(string observe_name)
Definition: observe.cc:305
unsigned int snap_dim_
Definition: observe.hh:106
const BIHTree & get_bih_tree()
Definition: mesh.cc:764
#define MPI_Comm_rank
Definition: mpi.h:236
Dedicated class for storing path to input and output files.
Definition: file_path.hh:48
void compute_field_values(Field< spacedim, Value > &field)
Definition: observe.cc:261
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:86
const double epsilon
Definition: mathfce.h:23
arma::vec3 global_coords_
Global coordinates of the observation point.
Definition: observe.hh:122
arma::vec local_coords_
Local (barycentric) coordinates on the element.
Definition: observe.hh:125
void update_projection(unsigned int i_elm, arma::vec local_coords, arma::vec3 global_coords)
Definition: observe.cc:82
double time_unit_seconds_
Time unit in seconds.
Definition: observe.hh:212
bool is_in_region_set(const RegionSet &set) const
Definition: region.cc:46
FieldCommon & name(const string &name)
Definition: field_common.hh:86
#define MPI_COMM_WORLD
Definition: mpi.h:123
string snap_region_name_
Definition: observe.hh:111
Record type proxy class.
Definition: type_record.hh:171
unsigned int element_idx_
Final element of the observe point. The index in the mesh.
Definition: observe.hh:119
void snap(Mesh &mesh)
Definition: observe.cc:121
#define INSTANCE_DIM(dim)
Definition: observe.cc:293
Class for declaration of the input data that are in string format.
Definition: type_base.hh:568
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:45
~Observe()
Destructor, must close the file.
Definition: observe.cc:255
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:219
void store_value(unsigned int idx, const Value &value)
Definition: output_data.cc:114
void output_time_frame(double time)
Definition: observe.cc:317
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:221