Flow123d  3.9.1-60c7e5c
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 "mesh/accessors.hh"
24 #include "io/observe.hh"
25 #include "io/element_data_cache.hh"
26 #include "fem/mapping_p1.hh"
27 #include "tools/time_governor.hh"
28 
29 
30 namespace IT = Input::Type;
31 
32 
33 /**
34  * Helper class allows to work with ObservePoint above elements with different dimensions
35  *
36  * Allows:
37  * - calculate projection of points by dimension
38  * - snap to subelement
39  */
40 template<unsigned int dim>
42 public:
43  /// Constructor
45 
46  ObservePointData projection(arma::vec3 input_point, unsigned int i_elm, ElementAccessor<3> elm) {
47  arma::mat::fixed<3, dim+1> elm_map = MappingP1<dim,3>::element_map(elm);
48  arma::vec::fixed<dim+1> projection = MappingP1<dim,3>::project_real_to_unit(input_point, elm_map);
50 
51  ObservePointData data;
52  data.element_idx_ = i_elm;
53  data.local_coords_ = projection.rows(1, elm.dim());
54  data.global_coords_ = elm_map*projection;
55  data.distance_ = arma::norm(data.global_coords_ - input_point, 2);
56  data.proc_ = elm.proc();
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 = MappingP1<dim,3>::element_map(elm);
81  observe_data.global_coords_ = elm_map * RefElement<dim>::local_to_bary(observe_data.local_coords_);
82  }
83 
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.\n"
109  "The actual observation element and the observation point on it is determined as follows:\n\n"
110  "1. Find an initial element containing the initial point. If no such element exists, we report an error.\n"
111  "2. Use BFS (Breadth-first search) starting from the inital element to find the 'observe element'. The observe element is the closest element.\n"
112  "3. Find the closest projection of the inital point on the observe element and snap this projection according to the ``snap_dim``.\n")
113  .allow_auto_conversion("point")
114  .declare_key("name", IT::String(),
116  "Default name have the form 'obs_<id>', where 'id' "
117  "is the rank of the point on the input."),
118  "Optional point name, which has to be unique.\n"
119  "Any string that is a valid YAML key in record without any quoting can be used, however, "
120  "using just alpha-numerical characters, and underscore instead of the space, is recommended."
121  )
123  "Initial point for the observe point search.")
124  .declare_key("snap_dim", IT::Integer(0, 4), IT::Default("4"),
125  "The dimension of the sub-element to which center we snap. For value 4 no snapping is done. "
126  "For values 0 up to 3 the element containing the initial point is found and then the observe"
127  "point is snapped to the nearest center of the sub-element of the given dimension. "
128  "E.g. for dimension 2 we snap to the nearest center of the face of the initial element."
129  )
130  .declare_key("snap_region", IT::String(), IT::Default("\"ALL\""),
131  "The region of the initial element for snapping. Without snapping we make a projection to the initial element.")
132  .declare_key("search_radius", IT::Double(0.0),
133  IT::Default::read_time("Maximal distance of the observe point from the mesh relative to the mesh diameter. "),
134  "Global value is defined in mesh record by the key global_snap_radius.")
135  .close();
136 }
137 
139 {}
140 
141 
142 ObservePoint::ObservePoint(Input::Record in_rec, Mesh &mesh, unsigned int point_idx)
143 {
144  in_rec_ = in_rec;
145 
146  string default_label = string("obs_") + std::to_string(point_idx);
147  name_ = in_rec.val<string>("name", default_label );
148 
149  vector<double> tmp_coords;
150  in_rec.val<Input::Array>("point").copy_to(tmp_coords);
151  input_point_= arma::vec(tmp_coords);
152 
153  snap_dim_ = in_rec.val<unsigned int>("snap_dim");
154 
155  snap_region_name_ = in_rec.val<string>("snap_region");
156 
157  const BoundingBox &main_box = mesh.get_bih_tree().tree_box();
158  double max_mesh_size = arma::max(main_box.max() - main_box.min());
159  max_search_radius_ = in_rec_.val<double>("search_radius", mesh.global_snap_radius()) * max_mesh_size;
160 }
161 
162 
163 
165  return observe_data_.distance_ < numeric_limits<double>::infinity();
166 }
167 
168 
169 
171 {
173  switch (elm.dim()) {
174  case 1:
175  {
178  break;
179  }
180  case 2:
181  {
184  break;
185  }
186  case 3:
187  {
190  break;
191  }
192  default: ASSERT_PERMANENT(false).error("Clipping supported only for dim=1,2,3.");
193  }
194 }
195 
196 
197 
200  if (region_set.size() == 0)
201  THROW( RegionDB::ExcUnknownSet() << RegionDB::EI_Label(snap_region_name_) << in_rec_.ei_address() );
202 
203 
204  const BIHTree &bih_tree=mesh.get_bih_tree();
205  vector<unsigned int> candidate_list;
206  std::unordered_set<unsigned int> closed_elements(1023);
207  std::priority_queue< ObservePointData, std::vector<ObservePointData>, CompareByDist > candidate_queue;
208 
209  // search for the initial element
210  auto projected_point = bih_tree.tree_box().project_point(input_point_);
211  bih_tree.find_point(projected_point, candidate_list, true);
212 
213  // closest element
214  ObservePointData min_observe_point_data;
215 
216  for (unsigned int i_candidate=0; i_candidate<candidate_list.size(); ++i_candidate) {
217  unsigned int i_elm=candidate_list[i_candidate];
218  ElementAccessor<3> elm = mesh.element_accessor(i_elm);
219 
220  // project point, add candidate to queue
221  auto observe_data = point_projection( i_elm, elm );
222 
223  // save the closest element for later diagnostic
224  if(observe_data.distance_ < min_observe_point_data.distance_)
225  min_observe_point_data = observe_data;
226 
227  // queue only the elements in the maximal search radius
228  if (observe_data.distance_ <= max_search_radius_)
229  candidate_queue.push(observe_data);
230  closed_elements.insert(i_elm);
231  }
232 
233  // no candidates found -> exception
234  if (candidate_queue.empty()) {
235  THROW(ExcNoObserveElementCandidates()
236  << EI_PointName(name_)
237  << EI_Point(input_point_)
238  << EI_ClosestEle(min_observe_point_data));
239  }
240 
241  while (!candidate_queue.empty())
242  {
243  auto candidate_data = candidate_queue.top();
244  candidate_queue.pop();
245 
246  unsigned int i_elm=candidate_data.element_idx_;
247  ElementAccessor<3> elm = mesh.element_accessor(i_elm);
248 
249  // test if candidate is in region and update projection
250  if (elm.region().is_in_region_set(region_set)) {
251  ASSERT_LE(candidate_data.distance_, observe_data_.distance_).error();
252 
253  observe_data_.distance_ = candidate_data.distance_;
254  observe_data_.element_idx_ = candidate_data.element_idx_;
255  observe_data_.local_coords_ = candidate_data.local_coords_;
256  observe_data_.global_coords_ = candidate_data.global_coords_;
257  observe_data_.proc_ = candidate_data.proc_;
258  break;
259  }
260 
261  // add candidates to queue
262  for (unsigned int n=0; n < elm->n_nodes(); n++)
263  for(unsigned int i_node_ele : mesh.node_elements()[elm.node(n).idx()]) {
264  if (closed_elements.find(i_node_ele) == closed_elements.end()) {
265  ElementAccessor<3> neighbor_elm = mesh.element_accessor(i_node_ele);
266  auto observe_data = point_projection( i_node_ele, neighbor_elm );
267  if (observe_data.distance_ <= max_search_radius_)
268  candidate_queue.push(observe_data);
269  closed_elements.insert(i_node_ele);
270  }
271  }
272  }
273 
274  if (! have_observe_element()) {
275  THROW(ExcNoObserveElement()
276  << EI_RegionName(snap_region_name_)
277  << EI_PointName(name_)
278  << EI_Point(input_point_)
279  << EI_ClosestEle(min_observe_point_data));
280  }
281  snap( mesh );
283  double dist = arma::norm(elm.centre() - input_point_, 2);
284  double elm_norm = arma::norm(elm.bounding_box().max() - elm.bounding_box().min(), 2);
285  if (dist > 2*elm_norm)
286  WarningOut().fmt("Observe point ({}) is too distant from the mesh.\n", name_);
287 }
288 
289 
290 
291 
292 void ObservePoint::output(ostream &out, unsigned int indent_spaces, unsigned int precision)
293 {
294  out << setw(indent_spaces) << "" << "- name: " << name_ << endl;
295  out << setw(indent_spaces) << "" << " init_point: " << field_value_to_yaml(input_point_, precision) << endl;
296  out << setw(indent_spaces) << "" << " snap_dim: " << snap_dim_ << endl;
297  out << setw(indent_spaces) << "" << " snap_region: " << snap_region_name_ << endl;
298  out << setw(indent_spaces) << "" << " observe_point: " << field_value_to_yaml(observe_data_.global_coords_, precision) << endl;
299 }
300 
301 
302 
304  switch (elm.dim()) {
305  case 1:
306  {
308  return ph.projection(input_point_, i_elm, elm);
309  break;
310  }
311  case 2:
312  {
314  return ph.projection(input_point_, i_elm, elm);
315  break;
316  }
317  case 3:
318  {
320  return ph.projection(input_point_, i_elm, elm);
321  break;
322  }
323  default:
324  ASSERT_PERMANENT(false).error("Invalid element dimension!");
325  }
326 
327  return ObservePointData(); // Should not happen.
328 }
329 
330 
331 
332 
333 /*******************************************************************
334  * implementation of Observe
335  */
336 
337 const unsigned int Observe::max_observe_value_time = 1000;
338 
339 
340 Observe::Observe(string observe_name, Mesh &mesh, Input::Array in_array,
341  unsigned int precision, const std::shared_ptr<TimeUnitConversion>& time_unit_conv)
342 : observe_name_(observe_name),
343  precision_(precision),
344  time_unit_conversion_(time_unit_conv),
345  point_ds_(nullptr),
346  observe_time_idx_(0)
347 {
349  observe_values_time_.push_back(numeric_limits<double>::signaling_NaN());
350 
351  unsigned int global_point_idx=0, local_point_idx=0;
352 
353  // in_rec is Output input record.
354  for(auto it = in_array.begin<Input::Record>(); it != in_array.end(); ++it) {
355  ObservePoint point(*it, mesh, points_.size());
356  point.find_observe_point(mesh);
357  point.observe_data_.global_idx_ = global_point_idx++;
358  if (point.observe_data_.proc_ == mesh.get_el_ds()->myp()) {
359  point.observe_data_.local_idx_ = local_point_idx++;
360  point_4_loc_.push_back(point.observe_data_.global_idx_);
361  }
362  else
363  point.observe_data_.local_idx_ = -1;
364  points_.push_back( point );
366  }
367  // make local to global map, distribution
368  point_ds_ = new Distribution(Observe::max_observe_value_time * point_4_loc_.size(), PETSC_COMM_WORLD);
369 
370  // make indices unique
371  std::sort(observed_element_indices_.begin(), observed_element_indices_.end());
372  auto last = std::unique(observed_element_indices_.begin(), observed_element_indices_.end());
374 
375  if (points_.size() == 0) return;
377  if (rank_==0) {
378  FilePath observe_file_path(observe_name_ + "_observe.yaml", FilePath::output_file);
379  try {
380  observe_file_path.open_stream(observe_file_);
381  //observe_file_.setf(std::ios::scientific);
382  observe_file_.precision(this->precision_);
383 
384  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, in_array)
385  output_header();
386  }
387 }
388 
390  if (points_.size()>0 && observe_field_values_.size()>0)
391  flush_values();
392 
393  observe_file_.close();
394  if (point_ds_!=nullptr) delete point_ds_;
395 }
396 
397 
398 template <typename T>
399 ElementDataCache<T> & Observe::prepare_compute_data(std::string field_name, double field_time, unsigned int n_rows,
400  unsigned int n_cols)
401 {
402  double time_unit_seconds = time_unit_conversion_->get_coef();
404  observe_values_time_[observe_time_idx_] = field_time / time_unit_seconds;
405  else
406  ASSERT(fabs(field_time / time_unit_seconds - observe_values_time_[observe_time_idx_]) < 2*numeric_limits<double>::epsilon())
408 
409  OutputDataFieldMap::iterator it=observe_field_values_.find(field_name);
410  if (it == observe_field_values_.end()) {
411  observe_field_values_[field_name]
412  = std::make_shared< ElementDataCache<T> >(field_name, n_rows * n_cols, point_ds_->lsize());
413  it=observe_field_values_.find(field_name);
414  }
415  return dynamic_cast<ElementDataCache<T> &>(*(it->second));
416 }
417 
418 // explicit instantiation of template method
419 #define OBSERVE_PREPARE_COMPUTE_DATA(TYPE) \
420 template ElementDataCache<TYPE> & Observe::prepare_compute_data<TYPE>(std::string field_name, double field_time, \
421  unsigned int n_rows, unsigned int n_cols)
422 
424 OBSERVE_PREPARE_COMPUTE_DATA(unsigned int);
426 
427 
429  unsigned int indent = 2;
430  observe_file_ << "# Observation file: " << observe_name_ << endl;
431  observe_file_ << "time_unit: " << time_unit_conversion_->get_unit_string() << endl;
432  observe_file_ << "time_unit_in_seconds: " << time_unit_conversion_->get_coef() << endl;
433  observe_file_ << "points:" << endl;
434  for(auto &point : points_)
435  point.output(observe_file_, indent, precision_);
436  observe_file_ << "data:" << endl;
437 
438 }
439 
442  for (unsigned int i=0; i<Observe::max_observe_value_time; ++i)
443  for (unsigned int j=0; j<point_4_loc_.size(); ++j) local_to_global[i*point_4_loc_.size()+j] = i*points_.size()+point_4_loc_[j];
444 
445  for(auto &field_data : observe_field_values_) {
446  auto serial_data = field_data.second->gather(point_ds_, &(local_to_global[0]));
447  if (rank_==0) field_data.second = serial_data;
448  }
449 
450  if (rank_ == 0) {
451  unsigned int indent = 2;
452  DebugOut() << "Observe::output_time_frame WRITE\n";
453  for (unsigned int i_time=0; i_time<observe_time_idx_; ++i_time) {
454  observe_file_ << setw(indent) << "" << "- time: " << observe_values_time_[i_time] << endl;
455  for(auto &field_data : observe_field_values_) {
456  observe_file_ << setw(indent) << "" << " " << field_data.second->field_input_name() << ": ";
457  field_data.second->print_yaml_subarray(observe_file_, precision_, i_time*points_.size(), (i_time+1)*points_.size());
458  observe_file_ << endl;
459  }
460  }
461  }
462 
463  observe_values_time_.clear();
465  observe_values_time_.push_back(numeric_limits<double>::signaling_NaN());
466  observe_time_idx_ = 0;
467 }
468 
469 void Observe::output_time_frame(bool flush) {
470  if (points_.size() == 0) return;
471 
472  if ( ! no_fields_warning ) {
473  no_fields_warning=true;
474  // check that observe fields are set
476  // first call and no fields
477  ASSERT(observe_field_values_.size() == 0);
478  WarningOut() << "No observe fields for the observation stream: " << observe_name_ << endl;
479  }
480  }
481 
483  ASSERT(observe_field_values_.size() == 0);
484  return;
485  }
486 
489  flush_values();
490  } else {
491  observe_values_time_.push_back( numeric_limits<double>::signaling_NaN() );
492  }
493 
494 }
495 
497  auto bgn_it = make_iter<ObservePointAccessor>( ObservePointAccessor(this, 0) );
498  auto end_it = make_iter<ObservePointAccessor>( ObservePointAccessor(this, point_4_loc_.size()) );
499  return Range<ObservePointAccessor>(bgn_it, end_it);
500 }
501 
ElementAccessor::dim
unsigned int dim() const
Definition: accessors.hh:188
Observe::ObservePointAccessor
friend class ObservePointAccessor
Definition: observe.hh:314
Observe::precision_
unsigned int precision_
Precision of float output.
Definition: observe.hh:298
element_data_cache.hh
Armor::vec
ArmaVec< double, N > vec
Definition: armor.hh:885
bih_tree.hh
Distribution::lsize
unsigned int lsize(int proc) const
get local size
Definition: distribution.hh:115
time_governor.hh
Basic time management class.
RefElement
Definition: ref_element.hh:339
MeshBase::region_db
const RegionDB & region_db() const
Definition: mesh.h:175
Distribution::myp
unsigned int myp() const
get my processor
Definition: distribution.hh:107
fmt::to_string
std::string to_string(const T &value)
Definition: string.h:29
ElementDataCache
Definition: element_data_cache.hh:44
ASSERT
#define ASSERT(expr)
Definition: asserts.hh:351
MeshBase::get_el_ds
Distribution * get_el_ds() const
Definition: mesh.h:119
ObservePointData::local_idx_
LongIdx local_idx_
Local index on actual process of the observe point.
Definition: observe.hh:66
Observe::point_ds_
Distribution * point_ds_
Parallel distribution of observe points.
Definition: observe.hh:309
Input::Type::Integer
Class for declaration of the integral input data.
Definition: type_base.hh:483
Input::Record::val
const Ret val(const string &key) const
Definition: accessors_impl.hh:31
ObservePointData::distance_
double distance_
Definition: observe.hh:57
CompareByDist::operator()
bool operator()(const ObservePointData &lhs, const ObservePointData &rhs) const
Definition: observe.cc:96
ProjectionHandler::ProjectionHandler
ProjectionHandler()
Constructor.
Definition: observe.cc:44
FilePath
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
Input::Type::Double
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
Observe::local_range
Range< ObservePointAccessor > local_range() const
Returns local range of observe points.
Definition: observe.cc:496
BIHTree
Class for O(log N) lookup for intersections with a set of bounding boxes.
Definition: bih_tree.hh:38
THROW
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
std::vector< double >
MappingP1::project_real_to_unit
static BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map)
Definition: mapping_p1.cc:69
Input::Type::Default::read_time
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
Observe::observe_values_time_
std::vector< double > observe_values_time_
Common evaluation time of the fields for single time frame.
Definition: observe.hh:289
ElementAccessor< 3 >
field_value_to_yaml
string field_value_to_yaml(const T &mat, unsigned int prec)
Definition: armadillo_tools.cc:138
arma::vec3
Definition: doxy_dummy_defs.hh:17
Observe::no_fields_warning
bool no_fields_warning
Definition: observe.hh:303
FilePath::open_stream
void open_stream(Stream &stream) const
Definition: file_path.cc:211
Input::Array::begin
Iterator< ValueType > begin() const
Definition: accessors_impl.hh:145
OBSERVE_PREPARE_COMPUTE_DATA
#define OBSERVE_PREPARE_COMPUTE_DATA(TYPE)
Definition: observe.cc:419
MPI_Comm_rank
#define MPI_Comm_rank
Definition: mpi.h:236
BoundingBox
Bounding box in 3d ambient space.
Definition: bounding_box.hh:53
ObservePoint::have_observe_element
bool have_observe_element()
Definition: observe.cc:164
Input::Type::Default
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Distribution
Definition: distribution.hh:50
ObservePointData
Definition: observe.hh:40
ProjectionHandler
Definition: observe.cc:41
Region::is_in_region_set
bool is_in_region_set(const RegionSet &set) const
Definition: region.cc:45
ObservePoint::get_input_type
static const Input::Type::Record & get_input_type()
Definition: observe.cc:107
ElementAccessor::bounding_box
BoundingBox bounding_box() const
Definition: accessors_impl.hh:142
armadillo_tools.hh
ObservePoint::snap_dim_
unsigned int snap_dim_
Definition: observe.hh:168
accessors.hh
BIHTree::find_point
void find_point(const Space< 3 >::Point &point, std::vector< unsigned int > &result_list, bool full_list=false) const
Definition: bih_tree.cc:287
Observe::flush_values
void flush_values()
Effectively writes the data into the observe stream.
Definition: observe.cc:440
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
Mesh::get_bih_tree
const BIHTree & get_bih_tree()
Getter for BIH. Creates and compute BIH at first call.
Definition: mesh.cc:1068
ASSERT_PERMANENT
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:348
BIHTree::tree_box
const BoundingBox & tree_box() const
Definition: bih_tree.cc:229
ObservePointData::proc_
unsigned int proc_
Actual process of the observe point.
Definition: observe.hh:60
ObservePoint::output
void output(ostream &out, unsigned int indent_spaces, unsigned int precision)
Definition: observe.cc:292
MeshBase::node_elements
const vector< vector< unsigned int > > & node_elements()
Definition: mesh.cc:1134
ObservePoint::snap
void snap(Mesh &mesh)
Definition: observe.cc:170
Input::Type::Record::allow_auto_conversion
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:133
accessors.hh
INPUT_CATCH
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:63
Observe::output_header
void output_header()
Definition: observe.cc:428
MeshBase::element_accessor
ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
Definition: mesh.cc:866
Input::Type::Default::obligatory
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
Observe::rank_
int rank_
Definition: observe.hh:277
Input::Record::ei_address
EI_Address ei_address() const
Definition: accessors.cc:178
ObservePoint
Definition: observe.hh:75
RegionDB::get_region_set
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:328
Observe::Observe
Observe(string observe_name, Mesh &mesh, Input::Array in_array, unsigned int precision, const std::shared_ptr< TimeUnitConversion > &time_unit_conv)
Definition: observe.cc:340
Input::Type::Record::declare_key
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:503
mesh.h
Observe::observe_file_
std::ofstream observe_file_
Output file stream.
Definition: observe.hh:295
fmt::internal::isnan
DummyInt isnan(...)
Definition: format.h:306
ElementAccessor::proc
unsigned int proc() const
Definition: accessors.hh:224
ElementAccessor::node
NodeAccessor< 3 > node(unsigned int ni) const
Definition: accessors.hh:230
ProjectionHandler::projection
ObservePointData projection(arma::vec3 input_point, unsigned int i_elm, ElementAccessor< 3 > elm)
Definition: observe.cc:46
Input::Type::Record::close
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
Input::Type
Definition: balance.hh:41
Observe::output_time_frame
void output_time_frame(bool flush)
Definition: observe.cc:469
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
ProjectionHandler::snap_to_subelement
void snap_to_subelement(ObservePointData &observe_data, ElementAccessor< 3 > elm, unsigned int snap_dim)
Definition: observe.cc:64
observe.hh
ObservePoint::snap_region_name_
string snap_region_name_
Definition: observe.hh:173
Observe::observed_element_indices_
std::vector< unsigned int > observed_element_indices_
Elements of the o_points.
Definition: observe.hh:282
MappingP1::element_map
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
input_type.hh
Observe::~Observe
~Observe()
Destructor, must close the file.
Definition: observe.cc:389
Mesh
Definition: mesh.h:362
MappingP1::clip_to_element
static BaryPoint clip_to_element(BaryPoint &barycentric)
Definition: mapping_p1.cc:91
Observe::prepare_compute_data
ElementDataCache< T > & prepare_compute_data(std::string field_name, double field_time, unsigned int n_rows, unsigned int n_cols)
Definition: observe.cc:399
Range
Range helper class.
Definition: range_wrapper.hh:65
ObservePointData::local_coords_
arma::vec local_coords_
Local (barycentric) coordinates on the element.
Definition: observe.hh:53
Input::Type::Array
Class for declaration of inputs sequences.
Definition: type_base.hh:339
Input::Type::String
Class for declaration of the input data that are in string format.
Definition: type_base.hh:582
Input::Array
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
global_defs.h
Global macros to enhance readability and debugging, general constants.
Observe::max_observe_value_time
static const unsigned int max_observe_value_time
Maximal size of observe values times vector.
Definition: observe.hh:274
WarningOut
#define WarningOut()
Macro defining 'warning' record of log.
Definition: logger.hh:278
ElementAccessor::region
Region region() const
Definition: accessors.hh:198
ASSERT_LE
#define ASSERT_LE(a, b)
Definition of comparative assert macro (Less or Equal) only for debug mode.
Definition: asserts.hh:309
Observe::observe_field_values_
OutputDataFieldMap observe_field_values_
Stored field values.
Definition: observe.hh:285
MPI_COMM_WORLD
#define MPI_COMM_WORLD
Definition: mpi.h:123
region.hh
ObservePoint::point_projection
ObservePointData point_projection(unsigned int i_elm, ElementAccessor< 3 > elm)
Project point to given element by dimension of this element.
Definition: observe.cc:303
ObservePoint::observe_data_
ObservePointData observe_data_
Helper object stored projection data.
Definition: observe.hh:184
Observe::point_4_loc_
std::vector< LongIdx > point_4_loc_
Index set assigning to local point index its global index.
Definition: observe.hh:306
Observe::observe_name_
std::string observe_name_
Definition: observe.hh:292
FilePath::output_file
@ output_file
Definition: file_path.hh:69
ObservePointData::element_idx_
unsigned int element_idx_
Final element of the observe point. The index in the mesh.
Definition: observe.hh:44
BoundingBox::max
const Point & max() const
Definition: bounding_box.hh:135
ObservePoint::max_search_radius_
double max_search_radius_
Definition: observe.hh:178
mapping_p1.hh
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
Observe::time_unit_conversion_
std::shared_ptr< TimeUnitConversion > time_unit_conversion_
Time unit conversion object.
Definition: observe.hh:300
RefElement::local_to_bary
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
Definition: ref_element.cc:187
Mesh::global_snap_radius
double global_snap_radius() const
Maximal distance of observe point from Mesh relative to its size.
Definition: mesh.cc:1077
Observe::observe_time_idx_
unsigned int observe_time_idx_
Index of actual (last) time in observe_values_time_ vector.
Definition: observe.hh:312
DebugOut
#define DebugOut()
Macro defining 'debug' record of log.
Definition: logger.hh:284
CompareByDist
Definition: observe.cc:94
Observe::points_
std::vector< ObservePoint > points_
Full information about observe points.
Definition: observe.hh:280
ObservePointData::global_idx_
LongIdx global_idx_
Global index of the observe point.
Definition: observe.hh:63
ObservePointData::global_coords_
arma::vec3 global_coords_
Global coordinates of the observation point.
Definition: observe.hh:50
ObservePoint::in_rec_
Input::Record in_rec_
Index in the input array.
Definition: observe.hh:159
BoundingBox::project_point
Point project_point(const Point &point) const
Definition: bounding_box.hh:286
Input::Array::end
IteratorBase end() const
Definition: accessors_impl.hh:157
BoundingBox::min
const Point & min() const
Definition: bounding_box.hh:128
ObservePoint::find_observe_point
void find_observe_point(Mesh &mesh)
Definition: observe.cc:198
ObservePoint::ObservePoint
ObservePoint()
Definition: observe.cc:138
ObservePoint::name_
std::string name_
Observation point name.
Definition: observe.hh:162
Element::n_nodes
unsigned int n_nodes() const
Definition: elements.h:125
ObservePoint::input_point_
arma::vec3 input_point_
Input coordinates of the initial position of the observation point.
Definition: observe.hh:181
ElementAccessor::centre
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
Definition: accessors_impl.hh:81