Flow123d  release_3.0.0-1159-g9e13a92
element_data_cache.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file element_data_cache.cc
15  * @brief
16  */
17 
18 
19 #include <limits>
20 #include <ostream>
21 #include "io/element_data_cache.hh"
22 #include "io/msh_basereader.hh"
23 #include "la/distribution.hh"
25 #include "system/system.hh"
26 #include "system/tokenizer.hh"
27 #include "boost/lexical_cast.hpp"
28 
29 
30 
31 template <typename T>
34  check_scale_data_(CheckScaleData::none) {}
35 
36 
37 template <typename T>
38 ElementDataCache<T>::ElementDataCache(std::string field_name, double time, unsigned int size_of_cache, unsigned int row_vec_size)
40 {
41  this->time_ = time;
42  this->field_input_name_ = field_name;
43  this->data_ = create_data_cache(size_of_cache, row_vec_size);
44 }
45 
46 
47 template <typename T>
48 ElementDataCache<T>::ElementDataCache(std::string field_name, unsigned int n_comp, unsigned int size)
50 {
51  this->set_vtk_type<T>();
52  this->field_name_ = field_name;
53  this->field_input_name_ = this->field_name_;
54 
55  this->n_values_ = size;
56  ASSERT_GT(n_comp, 0)(field_name).error("Output field returning variable size vectors. Try convert to MultiField.");
57  this->n_comp_ = n_comp;
58 
60 }
61 
62 
63 template <typename T>
65 
66 
67 template <typename T>
69  ASSERT_LT(component_idx, data_.size()).error("Index of component is out of range.\n");
70  return data_[component_idx];
71 }
72 
73 
74 template <typename T>
75 typename ElementDataCache<T>::CacheData ElementDataCache<T>::create_data_cache(unsigned int size_of_cache, unsigned int row_vec_size) {
76  typename ElementDataCache<T>::CacheData data_cache(size_of_cache);
77  for (unsigned int i=0; i<size_of_cache; ++i) {
78  typename ElementDataCache<T>::ComponentDataPtr row_vec = std::make_shared<std::vector<T>>();
79  row_vec->resize(row_vec_size, numeric_limits<T>::signaling_NaN());
80  data_cache[i] = row_vec;
81  }
82 
83  return data_cache;
84 }
85 
86 
87 template <typename T>
88 void ElementDataCache<T>::read_ascii_data(Tokenizer &tok, unsigned int n_components, unsigned int i_row) {
89  unsigned int idx;
90  for (unsigned int i_vec=0; i_vec<data_.size(); ++i_vec) {
91  idx = i_row * n_components;
92  std::vector<T> &vec = *( data_[i_vec].get() );
93  for (unsigned int i_col=0; i_col < n_components; ++i_col, ++idx) {
94  vec[idx] = boost::lexical_cast<T>(*tok);
95  ++tok;
96  }
97  }
98 }
99 
100 
101 template <typename T>
102 void ElementDataCache<T>::read_binary_data(std::istream &data_stream, unsigned int n_components, unsigned int i_row) {
103  unsigned int idx;
104  for (unsigned int i_vec=0; i_vec<data_.size(); ++i_vec) {
105  idx = i_row * n_components;
106  std::vector<T> &vec = *( data_[i_vec].get() );
107  for (unsigned int i_col=0; i_col < n_components; ++i_col, ++idx) {
108  data_stream.read(reinterpret_cast<char *>(&vec[idx]), sizeof(T));
109  }
110  }
111 }
112 
113 
114 /**
115  * Output data element on given index @p idx. Method for writing data
116  * to output stream.
117  *
118  * \note This method is used only by MSH file format.
119  */
120 template <typename T>
121 void ElementDataCache<T>::print_ascii(ostream &out_stream, unsigned int idx)
122 {
123  ASSERT_LT(idx, this->n_values_).error();
124  std::vector<T> &vec = *( this->data_[0].get() );
125  for(unsigned int i = n_comp_*idx; i < n_comp_*(idx+1); ++i )
126  out_stream << vec[i] << " ";
127 }
128 
129 /**
130  * \brief Print all data stored in output data
131  *
132  * TODO: indicate if the tensor data are output in column-first or raw-first order
133  * and possibly implement transposition. Set such property for individual file formats.
134  * Class OutputData stores always in raw-first order.
135  */
136 template <typename T>
137 void ElementDataCache<T>::print_ascii_all(ostream &out_stream)
138 {
139  std::vector<T> &vec = *( this->data_[0].get() );
140  for(unsigned int idx = 0; idx < this->n_values_; idx++) {
141  for(unsigned int i = n_comp_*idx; i < n_comp_*(idx+1); ++i )
142  out_stream << vec[i] << " ";
143  }
144 }
145 
146 
147 /// Prints the whole data vector into stream.
148 template <typename T>
149 void ElementDataCache<T>::print_binary_all(ostream &out_stream, bool print_data_size)
150 {
151  if (print_data_size) {
152  // write size of data
153  unsigned long long int data_byte_size = this->n_values_ * n_comp_ * sizeof(T);
154  out_stream.write(reinterpret_cast<const char*>(&data_byte_size), sizeof(unsigned long long int));
155  }
156  // write data
157  std::vector<T> &vec = *( this->data_[0].get() );
158  for(unsigned int idx = 0; idx < this->n_values_; idx++) {
159  for(unsigned int i = n_comp_*idx; i < n_comp_*(idx+1); ++i )
160  out_stream.write(reinterpret_cast<const char*>(&(vec[i])), sizeof(T));
161  }
162 }
163 
164 
165 template <typename T>
166 void ElementDataCache<T>::print_yaml_subarray(ostream &out_stream, unsigned int precision, unsigned int begin, unsigned int end)
167 {
168  out_stream << "[ ";
169  std::vector<T> &vec = *( this->data_[0].get() );
170  for(unsigned int idx = begin; idx < end; idx++) {
171  if (idx != begin) out_stream << " , ";
172  unsigned int vec_pos = n_comp_ * idx; // position of element value in data cache
173  switch (this->n_comp_) {
174  case NumCompValueType::N_SCALAR: {
175  out_stream << field_value_to_yaml( vec[vec_pos], precision );
176  break;
177  }
178  case NumCompValueType::N_VECTOR: {
179  typename arma::Col<T>::template fixed<3> vec_val;
180  for (unsigned int i=0; i<3; ++i, ++vec_pos)
181  vec_val(i) = vec[vec_pos];
182  out_stream << field_value_to_yaml( vec_val, precision );
183  break;
184  }
185  case NumCompValueType::N_TENSOR: {
186  typename arma::Mat<T>::template fixed<3,3> mat_val;
187  for (unsigned int i=0; i<3; ++i)
188  for (unsigned int j=0; j<3; ++j, ++vec_pos)
189  mat_val(i,j) = vec[vec_pos];
190  out_stream << field_value_to_yaml( mat_val, precision );
191  break;
192  }
193  }
194  }
195  out_stream << " ]";
196 }
197 
198 
199 template <typename T>
200 void ElementDataCache<T>::get_min_max_range(double &min, double &max)
201 {
202  min = std::numeric_limits<double>::max();
203  max = std::numeric_limits<double>::min();
204  std::vector<T> &vec = *( this->data_[0].get() );
205  for(unsigned int idx = 0; idx < this->n_values_; idx++) {
206  for(unsigned int i = n_comp_*idx; i < n_comp_*(idx+1); ++i ) {
207  if (vec[i] < min) min = vec[i];
208  if (vec[i] > max) max = vec[i];
209  }
210  }
211 }
212 
213 
214 /**
215  * Store data element of given data value under given index.
216  */
217 template <typename T>
218 void ElementDataCache<T>::store_value(unsigned int idx, const T * value) {
219  ASSERT_LT_DBG(idx, this->n_values_)(this->field_name_);
220  std::vector<T> &vec = *( this->data_[0].get() );
221  unsigned int vec_idx = idx*this->n_comp_;
222  for(unsigned int i = 0; i < this->n_comp_; i++, vec_idx++) {
223  vec[vec_idx] = value[i];
224  }
225 };
226 
227 /**
228  * Add value to given index
229  */
230 template <typename T>
231 void ElementDataCache<T>::add(unsigned int idx, const T * value) {
232  ASSERT_LT_DBG(idx, this->n_values_);
233  std::vector<T> &vec = *( this->data_[0].get() );
234  unsigned int vec_idx = idx*this->n_comp_;
235  for(unsigned int i = 0; i < this->n_comp_; i++, vec_idx++) {
236  vec[vec_idx] += value[i];
237  }
238 };
239 
240 /**
241  * Reset values at given index
242  */
243 template <typename T>
244 void ElementDataCache<T>::zero(unsigned int idx) {
245  ASSERT_LT_DBG(idx, this->n_values_);
246  std::vector<T> &vec = *( this->data_[0].get() );
247  unsigned int vec_idx = idx*this->n_comp_;
248  for(unsigned int i = 0; i < this->n_comp_; i++, vec_idx++) {
249  vec[vec_idx] = 0;
250  }
251 };
252 
253 /**
254  * Normalize values at given index
255  */
256 template <typename T>
257 void ElementDataCache<T>::normalize(unsigned int idx, unsigned int divisor) {
258  ASSERT_LT_DBG(idx, this->n_values_);
259  std::vector<T> &vec = *( this->data_[0].get() );
260  unsigned int vec_idx = idx*this->n_comp_;
261  for(unsigned int i = 0; i < this->n_comp_; i++, vec_idx++) {
262  vec[vec_idx] /= divisor;
263  }
264 };
265 
266 template <typename T>
267 CheckResult ElementDataCache<T>::check_values(double default_val, double lower_bound, double upper_bound) {
268  if (check_scale_data_ != CheckScaleData::none) return CheckResult::ok; // method is executed only once
270 
271  bool is_nan = false, out_of_limit = false;
272  for (unsigned int j=0; j<data_.size(); ++j) {
273  std::vector<T> &vec = *( this->data_[j].get() );
274  for(unsigned int i=0; i<vec.size(); ++i) {
275  if ( std::isnan(vec[i]) ) {
276  if ( std::isnan(default_val) ) is_nan = true;
277  else vec[i] = default_val;
278  }
279  if ( (vec[i] < lower_bound) || (vec[i] > upper_bound) ) out_of_limit = true;
280  }
281  }
282 
283  if (is_nan) return CheckResult::not_a_number;
284  else if (out_of_limit) return CheckResult::out_of_limits;
285  else return CheckResult::ok;
286 };
287 
288 template <typename T>
290  if (check_scale_data_ == CheckScaleData::scale) return; // method is executed only once
291  ASSERT_DBG(check_scale_data_ == CheckScaleData::check).warning("Data should be checked before scaling. Rather call 'check_values'!\n");
292 
293  for (unsigned int j=0; j<data_.size(); ++j) {
294  std::vector<T> &vec = *( this->data_[j].get() );
295  for(unsigned int i=0; i<vec.size(); ++i) {
296  vec[i] *= coef;
297  }
298  }
299 
300  check_scale_data_ = CheckScaleData::scale;
301 };
302 
303 
304 template <typename T>
305 std::shared_ptr< ElementDataCacheBase > ElementDataCache<T>::gather(Distribution *distr, LongIdx *local_to_global) {
306  std::shared_ptr< ElementDataCache<T> > gather_cache;
307  int rank = distr->myp();
308  int n_proc = distr->np();
309 
310  unsigned int n_global_data; // global number of data
311  int rec_starts[n_proc]; // displacement of first value that is received from each process
312  int rec_counts[n_proc]; // number of values that are received from each process
313  int *rec_indices_ids; // collective values of local to global indexes map of data
314  T *rec_data; // collective values of data
315 
316  // collects values of data vectors and local to global indexes map on each process
317  if (rank==0) {
318  for (int i=0; i<n_proc; ++i) {
319  rec_starts[i] = distr->begin(i);
320  rec_counts[i] = distr->lsize(i);
321  }
322  n_global_data = distr->size();
323  rec_indices_ids = new int [ n_global_data ];
324  }
325  MPI_Gatherv( local_to_global, distr->lsize(), MPI_INT, rec_indices_ids, rec_counts, rec_starts, MPI_INT, 0, MPI_COMM_WORLD);
326  if (rank==0) {
327  for (int i=0; i<n_proc; ++i) {
328  rec_starts[i] = this->n_comp()*rec_starts[i];
329  rec_counts[i] = this->n_comp()*rec_counts[i];
330  }
331  rec_data = new T [ this->n_comp() * n_global_data ];
332  }
333  auto &local_cache_vec = *( this->get_component_data(0).get() );
334  MPI_Gatherv( &local_cache_vec[0], this->n_comp()*distr->lsize(), this->mpi_data_type(), rec_data, rec_counts, rec_starts, this->mpi_data_type(), 0, MPI_COMM_WORLD) ;
335 
336  // create and fill serial cache
337  if (rank==0) {
338  gather_cache = std::make_shared<ElementDataCache<T>>(this->field_input_name_, (unsigned int)this->n_comp(), n_global_data);
339  auto &gather_vec = *( gather_cache->get_component_data(0).get() );
340  unsigned int i_global_coord; // counter over serial_mesh->nodes_ cache
341  for (unsigned int i=0; i<n_global_data; ++i) {
342  i_global_coord = this->n_comp() * rec_indices_ids[i];
343  for (unsigned int j=0; j<this->n_comp(); ++j) { //loop over coords
344  ASSERT_LT(i_global_coord+j, gather_vec.size());
345  gather_vec[ i_global_coord+j ] = rec_data[ this->n_comp()*i+j ];
346  }
347  }
348 
349  delete[] rec_indices_ids;
350  delete[] rec_data;
351  }
352 
353  return gather_cache;
354 }
355 
356 
357 template <typename T>
358 std::shared_ptr< ElementDataCacheBase > ElementDataCache<T>::element_node_cache_fixed_size(std::vector<unsigned int> &offset_vec) {
359  unsigned int n_elem = offset_vec.size();
360  std::shared_ptr< ElementDataCache<T> > elem_node_cache = std::make_shared<ElementDataCache<T>>(this->field_input_name_, 4*this->n_comp(), n_elem);
361  auto &data_out_vec = *( elem_node_cache->get_component_data(0).get() );
362  std::fill( data_out_vec.begin(), data_out_vec.end(), (T)0 );
363  auto &data_in_vec = *( this->get_component_data(0).get() );
364 
365  unsigned int i_node, i_old, i_new;
366  for (unsigned int i_el=0, i_conn=0; i_el<offset_vec.size(); i_el++) {
367  for(i_node=4*i_el; i_conn<offset_vec[i_el]; i_conn++, i_node++) {
368  i_old = i_conn*this->n_comp_;
369  i_new = i_node*this->n_comp_;
370  for(unsigned int i = 0; i < this->n_comp_; i++) {
371  ASSERT_LT(i_new+i, data_out_vec.size());
372  ASSERT_LT(i_old+i, data_in_vec.size());
373  data_out_vec[i_new+i] = data_in_vec[i_old+i];
374  }
375  }
376  }
377 
378  return elem_node_cache;
379 }
380 
381 
382 template <typename T>
383 std::shared_ptr< ElementDataCacheBase > ElementDataCache<T>::element_node_cache_optimize_size(std::vector<unsigned int> &offset_vec) {
384  std::shared_ptr< ElementDataCache<T> > elem_node_cache = std::make_shared<ElementDataCache<T>>(this->field_input_name_,
385  this->n_comp()/4, offset_vec[offset_vec.size()-1]);
386  auto &data_out_vec = *( elem_node_cache->get_component_data(0).get() );
387  auto &data_in_vec = *( this->get_component_data(0).get() );
388 
389  unsigned int i_node, i_old, i_new;
390  for (unsigned int i_el=0, i_conn=0; i_el<offset_vec.size(); i_el++) {
391  for(i_node=4*i_el; i_conn<offset_vec[i_el]; i_conn++, i_node++) {
392  i_old = i_node*elem_node_cache->n_comp_;
393  i_new = i_conn*elem_node_cache->n_comp_;
394  for(unsigned int i = 0; i < elem_node_cache->n_comp_; i++) {
395  ASSERT_LT(i_new+i, data_out_vec.size());
396  ASSERT_LT(i_old+i, data_in_vec.size());
397  data_out_vec[i_new+i] = data_in_vec[i_old+i];
398  }
399  }
400  }
401  return elem_node_cache;
402 }
403 
404 
405 template <typename T>
406 std::shared_ptr< ElementDataCacheBase > ElementDataCache<T>::compute_node_data(std::vector<unsigned int> &conn_vec, unsigned int data_size) {
407  ASSERT_EQ(conn_vec.size(), this->n_values());
408  unsigned int idx;
409 
410  // set output data to zero
411  std::shared_ptr< ElementDataCache<T> > node_cache = std::make_shared<ElementDataCache<T>>(this->field_input_name_, this->n_comp(), data_size);
412  std::vector<unsigned int> count(data_size, 0);
413  for (idx=0; idx < node_cache->n_values(); idx++)
414  node_cache->zero(idx);
415 
416  auto &data_in_vec = *( this->get_component_data(0).get() );
417  for (idx=0; idx < conn_vec.size(); idx++) {
418  ASSERT_LT(conn_vec[idx], node_cache->n_values());
419  ASSERT_LT(this->n_comp()*idx, data_in_vec.size());
420  node_cache->add( conn_vec[idx], &(data_in_vec[this->n_comp() * idx]) );
421  count[ conn_vec[idx] ]++;
422  }
423 
424  // Compute mean values at nodes
425  for(idx=0; idx < node_cache->n_values(); idx++)
426  node_cache->normalize(idx, count[idx]);
427 
428  return node_cache;
429 }
430 
431 
432 template<>
434  return MPI_DOUBLE;
435 };
436 
437 template<>
439  return MPI_INT;
440 };
441 
442 template<>
444  return MPI_UNSIGNED;
445 };
446 
447 
448 /// Access i-th element in the data vector.
449 template <class T>
451 {
452  std::vector<T> &vec = *( this->data_[0].get() );
453  ASSERT_DBG(i < vec.size());
454  return vec[i];
455 }
456 
457 
458 
459 // explicit instantiation of template class
460 template class ElementDataCache<unsigned int>;
461 template class ElementDataCache<int>;
462 template class ElementDataCache<double>;
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
unsigned int size() const
get global size
double time_
time step stored in cache
void read_binary_data(std::istream &data_stream, unsigned int n_components, unsigned int i_row) override
Implements ElementDataCacheBase::read_binary_data.
std::shared_ptr< std::vector< T > > ComponentDataPtr
CheckScaleData
Allow to hold sign, if data in cache is checked and scale (both can be executed only once) ...
#define MPI_Datatype
Definition: mpi.h:154
void print_ascii(ostream &out_stream, unsigned int idx) override
string field_value_to_yaml(const T &mat, unsigned int prec)
CheckScaleData check_scale_data_
Sign, if data in cache is checked and scale.
Some value(s) is set to NaN.
unsigned int n_values() const
std::shared_ptr< ElementDataCacheBase > compute_node_data(std::vector< unsigned int > &conn_vec, unsigned int data_size) override
Implements ElementDataCacheBase::compute_node_data.
void store_value(unsigned int idx, const T *value)
#define ASSERT_GT(a, b)
Definition of comparative assert macro (Greater Than)
Definition: asserts.hh:311
#define MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm)
Definition: mpi.h:549
T & operator[](unsigned int i)
Access i-th element in the data vector of 0th component.
Mat< double, N, 1 > vec
Definition: armor.hh:211
MPI_Datatype mpi_data_type()
Return MPI data type corresponding with template parameter of cache. Needs template specialization...
void read_ascii_data(Tokenizer &tok, unsigned int n_components, unsigned int i_row) override
Implements ElementDataCacheBase::read_ascii_data.
Some value(s) is out of limits.
static constexpr bool value
Definition: json.hpp:87
ElementDataCache()
Default constructor.
void get_min_max_range(double &min, double &max) override
void add(unsigned int idx, const T *value)
void print_yaml_subarray(ostream &out_stream, unsigned int precision, unsigned int begin, unsigned int end) override
DummyInt isnan(...)
Definition: format.h:306
std::shared_ptr< ElementDataCacheBase > gather(Distribution *distr, LongIdx *local_to_global) override
Implements ElementDataCacheBase::gather.
virtual ~ElementDataCache() override
Destructor of ElementDataCache.
static CacheData create_data_cache(unsigned int size_of_cache, unsigned int row_vec_size)
unsigned int begin(int proc) const
get starting local index
void zero(unsigned int idx)
void scale_data(double coef)
unsigned int np() const
get num of processors
Data is neither checked nor scaled.
CheckResult
Return type of method that checked data stored in ElementDataCache (NaN values, limits) ...
CheckResult check_values(double default_val, double lower_bound, double upper_bound)
#define MPI_DOUBLE
Definition: mpi.h:156
void normalize(unsigned int idx, unsigned int divisor)
unsigned int myp() const
get my processor
Support classes for parallel programing.
void print_binary_all(ostream &out_stream, bool print_data_size=true) override
Print all data stored in output data to appended binary format.
std::string field_input_name_
name of field stored in cache
std::shared_ptr< ElementDataCacheBase > element_node_cache_fixed_size(std::vector< unsigned int > &offset_vec) override
Implements ElementDataCacheBase::element_node_cache_fixed_size.
#define MPI_INT
Definition: mpi.h:160
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than)
Definition: asserts.hh:295
std::shared_ptr< ElementDataCacheBase > element_node_cache_optimize_size(std::vector< unsigned int > &offset_vec) override
Implements ElementDataCacheBase::element_node_cache_optimize_size.
ComponentDataPtr get_component_data(unsigned int component_idx)
Return vector of element data for get component.
#define ASSERT_DBG(expr)
Definition: asserts.hh:349
#define MPI_COMM_WORLD
Definition: mpi.h:123
void print_ascii_all(ostream &out_stream) override
Print all data stored in output data ro ascii format.
All values are not NaN and are in limits.
#define MPI_UNSIGNED
Definition: mpi.h:165
unsigned int n_comp() const
T check(T value)
Definition: format.h:312
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:327
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:299
unsigned int lsize(int proc) const
get local size