Flow123d  release_3.0.0-922-g9fcbba1
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  data.proc_ = elm.proc();
58 
59  return data;
60  }
61 
62  /**
63  * Snap local coords to the subelement. Called by the ObservePoint::snap method.
64  */
65  void snap_to_subelement(ObservePointData & observe_data, ElementAccessor<3> elm, unsigned int snap_dim)
66  {
67  if (snap_dim <= dim) {
68  double min_dist = 2.0; // on the ref element the max distance should be about 1.0, smaler then 2.0
69  arma::vec min_center;
70  for(auto &center : RefElement<dim>::centers_of_subelements(snap_dim))
71  {
72  double dist = arma::norm(center - observe_data.local_coords_, 2);
73  if ( dist < min_dist) {
74  min_dist = dist;
75  min_center = center;
76  }
77  }
78  observe_data.local_coords_ = min_center;
79  }
80 
81  arma::mat::fixed<3, dim+1> elm_map = mapping_.element_map(elm);
82  observe_data.global_coords_ = elm_map * RefElement<dim>::local_to_bary(observe_data.local_coords_);
83  }
84 
85 private:
86  /// Mapping object.
88 };
89 
90 template class ProjectionHandler<1>;
91 template class ProjectionHandler<2>;
92 template class ProjectionHandler<3>;
93 
94 
95 /**
96  * Helper struct, used as comparator of priority queue in ObservePoint::find_observe_point.
97  */
99 {
100  bool operator()(const ObservePointData& lhs, const ObservePointData& rhs) const
101  {
102  return lhs.distance_ > rhs.distance_;
103  }
104 };
105 
106 
107 /*******************************************************************
108  * implementation of ObservePoint
109  */
110 
112  return IT::Record("ObservePoint", "Specification of the observation point.\n"
113  "The actual observation element and the observation point on it is determined as follows:\n\n"
114  "1. Find an initial element containing the initial point. If no such element exists, we report an error.\n"
115  "2. Use BFS (Breadth-first search) starting from the inital element to find the 'observe element'. The observe element is the closest element.\n"
116  "3. Find the closest projection of the inital point on the observe element and snap this projection according to the ``snap_dim``.\n")
117  .allow_auto_conversion("point")
118  .declare_key("name", IT::String(),
120  "Default name have the form 'obs_<id>', where 'id' "
121  "is the rank of the point on the input."),
122  "Optional point name, which has to be unique.\n"
123  "Any string that is a valid YAML key in record without any quoting can be used, however, "
124  "using just alpha-numerical characters, and underscore instead of the space, is recommended."
125  )
127  "Initial point for the observe point search.")
128  .declare_key("snap_dim", IT::Integer(0, 4), IT::Default("4"),
129  "The dimension of the sub-element to which center we snap. For value 4 no snapping is done. "
130  "For values 0 up to 3 the element containing the initial point is found and then the observe"
131  "point is snapped to the nearest center of the sub-element of the given dimension. "
132  "E.g. for dimension 2 we snap to the nearest center of the face of the initial element."
133  )
134  .declare_key("snap_region", IT::String(), IT::Default("\"ALL\""),
135  "The region of the initial element for snapping. Without snapping we make a projection to the initial element.")
136  .declare_key("search_radius", IT::Double(0.0),
137  IT::Default::read_time("Maximal distance of observe point from Mesh relative to its size (bounding box). "),
138  "Global value is define in Mesh by the key global_snap_radius.")
139  .close();
140 }
141 
143 {}
144 
145 
146 ObservePoint::ObservePoint(Input::Record in_rec, Mesh &mesh, unsigned int point_idx)
147 {
148  in_rec_ = in_rec;
149 
150  string default_label = string("obs_") + std::to_string(point_idx);
151  name_ = in_rec.val<string>("name", default_label );
152 
153  vector<double> tmp_coords;
154  in_rec.val<Input::Array>("point").copy_to(tmp_coords);
155  input_point_= arma::vec(tmp_coords);
156 
157  snap_dim_ = in_rec.val<unsigned int>("snap_dim");
158 
159  snap_region_name_ = in_rec.val<string>("snap_region");
160 
161  const BoundingBox &main_box = mesh.get_bih_tree().tree_box();
162  double max_mesh_size = arma::max(main_box.max() - main_box.min());
163  max_search_radius_ = in_rec_.val<double>("search_radius", mesh.global_snap_radius()) * max_mesh_size;
164 }
165 
166 
167 
169  return observe_data_.distance_ < numeric_limits<double>::infinity();
170 }
171 
172 
173 
175 {
176  ElementAccessor<3> elm = mesh.element_accessor(observe_data_.element_idx_);
177  switch (elm.dim()) {
178  case 1:
179  {
181  ph.snap_to_subelement(observe_data_, elm, snap_dim_);
182  break;
183  }
184  case 2:
185  {
187  ph.snap_to_subelement(observe_data_, elm, snap_dim_);
188  break;
189  }
190  case 3:
191  {
193  ph.snap_to_subelement(observe_data_, elm, snap_dim_);
194  break;
195  }
196  default: ASSERT(false).error("Clipping supported only for dim=1,2,3.");
197  }
198 }
199 
200 
201 
203  RegionSet region_set = mesh.region_db().get_region_set(snap_region_name_);
204  if (region_set.size() == 0)
205  THROW( RegionDB::ExcUnknownSet() << RegionDB::EI_Label(snap_region_name_) << in_rec_.ei_address() );
206 
207 
208  const BIHTree &bih_tree=mesh.get_bih_tree();
209  vector<unsigned int> candidate_list;
210  std::unordered_set<unsigned int> closed_elements(1023);
211  std::priority_queue< ObservePointData, std::vector<ObservePointData>, CompareByDist > candidate_queue;
212 
213  // search for the initial element
214  auto projected_point = bih_tree.tree_box().project_point(input_point_);
215  bih_tree.find_point(projected_point, candidate_list, true);
216 
217  for (unsigned int i_candidate=0; i_candidate<candidate_list.size(); ++i_candidate) {
218  unsigned int i_elm=candidate_list[i_candidate];
219  ElementAccessor<3> elm = mesh.element_accessor(i_elm);
220 
221  // project point, add candidate to queue
222  auto observe_data = point_projection( i_elm, elm );
223  if (observe_data.distance_ <= max_search_radius_)
224  candidate_queue.push(observe_data);
225  closed_elements.insert(i_elm);
226  }
227 
228  while (!candidate_queue.empty())
229  {
230  auto candidate_data = candidate_queue.top();
231  candidate_queue.pop();
232 
233  unsigned int i_elm=candidate_data.element_idx_;
234  ElementAccessor<3> elm = mesh.element_accessor(i_elm);
235 
236  // test if candidate is in region and update projection
237  if (elm.region().is_in_region_set(region_set)) {
238  ASSERT_LE(candidate_data.distance_, observe_data_.distance_).error();
239 
240  observe_data_.distance_ = candidate_data.distance_;
241  observe_data_.element_idx_ = candidate_data.element_idx_;
242  observe_data_.local_coords_ = candidate_data.local_coords_;
243  observe_data_.global_coords_ = candidate_data.global_coords_;
244  observe_data_.proc_ = candidate_data.proc_;
245  break;
246  }
247 
248  // add candidates to queue
249  for (unsigned int n=0; n < elm->n_nodes(); n++)
250  for(unsigned int i_node_ele : mesh.node_elements()[elm.node_accessor(n).idx()]) {
251  if (closed_elements.find(i_node_ele) == closed_elements.end()) {
252  ElementAccessor<3> neighbor_elm = mesh.element_accessor(i_node_ele);
253  auto observe_data = point_projection( i_node_ele, neighbor_elm );
254  if (observe_data.distance_ <= max_search_radius_)
255  candidate_queue.push(observe_data);
256  closed_elements.insert(i_node_ele);
257  }
258  }
259  }
260 
261  if (! have_observe_element()) {
262  THROW(ExcNoObserveElement() << EI_RegionName(snap_region_name_) );
263  }
264  snap( mesh );
265  ElementAccessor<3> elm = mesh.element_accessor(observe_data_.element_idx_);
266  double dist = arma::norm(elm.centre() - input_point_, 2);
267  double elm_norm = arma::norm(elm.bounding_box().max() - elm.bounding_box().min(), 2);
268  if (dist > 2*elm_norm)
269  WarningOut().fmt("Observe point ({}) is too distant from the mesh.\n", name_);
270 }
271 
272 
273 
274 
275 void ObservePoint::output(ostream &out, unsigned int indent_spaces, unsigned int precision)
276 {
277  out << setw(indent_spaces) << "" << "- name: " << name_ << endl;
278  out << setw(indent_spaces) << "" << " init_point: " << field_value_to_yaml(input_point_, precision) << endl;
279  out << setw(indent_spaces) << "" << " snap_dim: " << snap_dim_ << endl;
280  out << setw(indent_spaces) << "" << " snap_region: " << snap_region_name_ << endl;
281  out << setw(indent_spaces) << "" << " observe_point: " << field_value_to_yaml(observe_data_.global_coords_, precision) << endl;
282 }
283 
284 
285 
287  switch (elm.dim()) {
288  case 1:
289  {
291  return ph.projection(input_point_, i_elm, elm);
292  break;
293  }
294  case 2:
295  {
297  return ph.projection(input_point_, i_elm, elm);
298  break;
299  }
300  case 3:
301  {
303  return ph.projection(input_point_, i_elm, elm);
304  break;
305  }
306  default:
307  ASSERT(false).error("Invalid element dimension!");
308  }
309 
310  return ObservePointData(); // Should not happen.
311 }
312 
313 
314 
315 
316 /*******************************************************************
317  * implementation of Observe
318  */
319 
320 const unsigned int Observe::max_observe_value_time = 1000;
321 
322 
323 Observe::Observe(string observe_name, Mesh &mesh, Input::Array in_array, unsigned int precision, std::string unit_str)
324 : observe_name_(observe_name),
325  precision_(precision),
326  point_ds_(nullptr),
327  observe_time_idx_(0)
328 {
330  observe_values_time_.push_back(numeric_limits<double>::signaling_NaN());
331 
332  unsigned int global_point_idx=0, local_point_idx=0;
333 
334  // in_rec is Output input record.
335  for(auto it = in_array.begin<Input::Record>(); it != in_array.end(); ++it) {
336  ObservePoint point(*it, mesh, points_.size());
337  point.find_observe_point(mesh);
338  point.observe_data_.global_idx_ = global_point_idx++;
339  if (point.observe_data_.proc_ == mesh.get_el_ds()->myp()) {
340  point.observe_data_.local_idx_ = local_point_idx++;
341  point_4_loc_.push_back(point.observe_data_.global_idx_);
342  }
343  else
344  point.observe_data_.local_idx_ = -1;
345  points_.push_back( point );
346  observed_element_indices_.push_back(point.observe_data_.element_idx_);
347  }
348  // make local to global map, distribution
349  point_ds_ = new Distribution(Observe::max_observe_value_time * point_4_loc_.size(), PETSC_COMM_WORLD);
350 
351  // make indices unique
352  std::sort(observed_element_indices_.begin(), observed_element_indices_.end());
353  auto last = std::unique(observed_element_indices_.begin(), observed_element_indices_.end());
355 
356  time_unit_str_ = unit_str;
358 
359  if (points_.size() == 0) return;
361  if (rank_==0) {
362  FilePath observe_file_path(observe_name_ + "_observe.yaml", FilePath::output_file);
363  try {
364  observe_file_path.open_stream(observe_file_);
365  //observe_file_.setf(std::ios::scientific);
366  observe_file_.precision(this->precision_);
367 
368  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, in_array)
369  output_header();
370  }
371 }
372 
374  observe_file_.close();
375  if (point_ds_!=nullptr) delete point_ds_;
376 }
377 
378 
379 template <typename T>
380 ElementDataCache<T> & Observe::prepare_compute_data(std::string field_name, double field_time, unsigned int n_rows,
381  unsigned int n_cols)
382 {
385  else
386  ASSERT(fabs(field_time / time_unit_seconds_ - observe_values_time_[observe_time_idx_]) < 2*numeric_limits<double>::epsilon())
387  (field_time)(observe_values_time_[observe_time_idx_]);
388 
389  OutputDataFieldMap::iterator it=observe_field_values_.find(field_name);
390  if (it == observe_field_values_.end()) {
391  observe_field_values_[field_name]
392  = std::make_shared< ElementDataCache<T> >(field_name, n_rows * n_cols, point_ds_->lsize());
393  it=observe_field_values_.find(field_name);
394  }
395  return dynamic_cast<ElementDataCache<T> &>(*(it->second));
396 }
397 
398 // explicit instantiation of template method
399 #define OBSERVE_PREPARE_COMPUTE_DATA(TYPE) \
400 template ElementDataCache<TYPE> & Observe::prepare_compute_data<TYPE>(std::string field_name, double field_time, \
401  unsigned int n_rows, unsigned int n_cols)
402 
404 OBSERVE_PREPARE_COMPUTE_DATA(unsigned int);
406 
407 
409  unsigned int indent = 2;
410  observe_file_ << "# Observation file: " << observe_name_ << endl;
411  observe_file_ << "time_unit: " << time_unit_str_ << endl;
412  observe_file_ << "time_unit_in_seconds: " << time_unit_seconds_ << endl;
413  observe_file_ << "points:" << endl;
414  for(auto &point : points_)
415  point.output(observe_file_, indent, precision_);
416  observe_file_ << "data:" << endl;
417 
418 }
419 
420 void Observe::output_time_frame(double time, bool flush) {
421  if (points_.size() == 0) return;
422 
423  if ( ! no_fields_warning ) {
424  no_fields_warning=true;
425  // check that observe fields are set
427  // first call and no fields
428  ASSERT(observe_field_values_.size() == 0);
429  WarningOut() << "No observe fields for the observation stream: " << observe_name_ << endl;
430  }
431  }
432 
434  ASSERT(observe_field_values_.size() == 0);
435  return;
436  }
437 
438  observe_time_idx_++;
439  if ( (observe_time_idx_==Observe::max_observe_value_time) || flush ) {
441  for (unsigned int i=0; i<Observe::max_observe_value_time; ++i)
442  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];
443 
444  for(auto &field_data : observe_field_values_) {
445  auto serial_data = field_data.second->gather(point_ds_, &(local_to_global[0]));
446  if (rank_==0) field_data.second = serial_data;
447  }
448 
449  if (rank_ == 0) {
450  unsigned int indent = 2;
451  for (unsigned int i_time=0; i_time<observe_time_idx_; ++i_time) {
452  observe_file_ << setw(indent) << "" << "- time: " << observe_values_time_[i_time] << endl;
453  for(auto &field_data : observe_field_values_) {
454  observe_file_ << setw(indent) << "" << " " << field_data.second->field_input_name() << ": ";
455  field_data.second->print_yaml_subarray(observe_file_, precision_, i_time*points_.size(), (i_time+1)*points_.size());
456  observe_file_ << endl;
457  }
458  }
459  }
460 
461  observe_values_time_.clear();
463  observe_values_time_.push_back(numeric_limits<double>::signaling_NaN());
464  observe_time_idx_ = 0;
465  } else {
466  observe_values_time_.push_back( numeric_limits<double>::signaling_NaN() );
467  }
468 
469 }
470 
472  auto bgn_it = make_iter<ObservePointAccessor>( ObservePointAccessor(this, 0) );
473  auto end_it = make_iter<ObservePointAccessor>( ObservePointAccessor(this, point_4_loc_.size()) );
474  return Range<ObservePointAccessor>(bgn_it, end_it);
475 }
476 
ObservePointData point_projection(unsigned int i_elm, ElementAccessor< 3 > elm)
Project point to given element by dimension of this element.
Definition: observe.cc:286
std::vector< ObservePoint > points_
Full information about observe points.
Definition: observe.hh:260
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:399
OutputDataFieldMap observe_field_values_
Stored field values.
Definition: observe.hh:265
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:323
vector< vector< unsigned int > > const & node_elements()
Definition: mesh.cc:999
std::ofstream observe_file_
Output file stream.
Definition: observe.hh:275
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
unsigned int proc_
Actual process of the observe point.
Definition: observe.hh:58
BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map) const
Definition: mapping_p1.cc:275
std::vector< LongIdx > point_4_loc_
Index set assigning to local point index its global index.
Definition: observe.hh:288
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:111
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:64
Range helper class.
Definition: mesh.h:52
void output_time_frame(double time, bool flush)
Definition: observe.cc:420
bool have_observe_element()
Definition: observe.cc:168
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:941
Definition: mesh.h:80
void output(ostream &out, unsigned int indent_spaces, unsigned int precision)
Definition: observe.cc:275
static const unsigned int max_observe_value_time
Maximal size of observe values times vector.
Definition: observe.hh:254
#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:42
std::string time_unit_str_
String representation of the time unit.
Definition: observe.hh:278
int rank_
Definition: observe.hh:257
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:272
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:51
unsigned int proc() const
Definition: accessors.hh:125
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:727
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:48
unsigned int precision_
Precision of float output.
Definition: observe.hh:282
void find_observe_point(Mesh &mesh)
Definition: observe.cc:202
DummyInt isnan(...)
Definition: format.h:306
bool operator()(const ObservePointData &lhs, const ObservePointData &rhs) const
Definition: observe.cc:100
const BoundingBox & tree_box() const
Definition: bih_tree.cc:229
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:408
std::vector< unsigned int > observed_element_indices_
Elements of the o_points.
Definition: observe.hh:262
const Point & max() const
double distance_
Definition: observe.hh:55
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:65
Range< ObservePointAccessor > local_range() const
Returns local range of observe points.
Definition: observe.cc:471
ElementDataCache< T > & prepare_compute_data(std::string field_name, double field_time, unsigned int n_rows, unsigned int n_cols)
Definition: observe.cc:380
Distribution * get_el_ds() const
Definition: mesh.h:168
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:932
#define MPI_Comm_rank
Definition: mpi.h:236
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
unsigned int myp() const
get my processor
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
friend class ObservePointAccessor
Definition: observe.hh:296
void find_point(const Space< 3 >::Point &point, std::vector< unsigned int > &result_list, bool full_list=false) const
Definition: bih_tree.cc:287
double time_unit_seconds_
Time unit in seconds.
Definition: observe.hh:280
std::vector< double > observe_values_time_
Common evaluation time of the fields for single time frame.
Definition: observe.hh:269
bool is_in_region_set(const RegionSet &set) const
Definition: region.cc:46
MappingP1< dim, 3 > mapping_
Mapping object.
Definition: observe.cc:87
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
unsigned int observe_time_idx_
Index of actual (last) time in observe_values_time_ vector.
Definition: observe.hh:294
const Point & min() const
Record type proxy class.
Definition: type_record.hh:182
Distribution * point_ds_
Parallel distribution of observe points.
Definition: observe.hh:291
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:174
bool no_fields_warning
Definition: observe.hh:285
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:373
unsigned int lsize(int proc) const
get local size