Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_Sparse3TensorPartition.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef STOKHOS_SPARSE_3_TENSOR_PARTITION_HPP
11 #define STOKHOS_SPARSE_3_TENSOR_PARTITION_HPP
12 
14 
15 #include "Teuchos_ArrayRCP.hpp"
16 #include "Teuchos_Array.hpp"
17 
18 namespace Stokhos {
19 
20  template <typename TupleType>
21  class RCB {
22  public:
25  typedef typename TupleType::id_type id_type;
26 
27  struct CoordCompare {
29  CoordCompare(const size_type& d_) : d(d_) {}
30  bool operator() (const TupleType& a, const TupleType& b) const {
31  return a(d) < b(d);
32  }
33  };
34 
35  struct Box {
41 
42  Box(const Teuchos::ArrayView<TupleType>& c) : coords(c.begin(), c.end()) {
45  }
46 
47  // Compute bounding box around points
49  xmin = coords[0](0); xmax = coords[0](0);
50  ymin = coords[0](1); ymax = coords[0](1);
51  zmin = coords[0](2); zmax = coords[0](2);
52  for (size_type i=0; i<coords.size(); ++i) {
53  coord_type x = coords[i](0);
54  coord_type y = coords[i](1);
55  coord_type z = coords[i](2);
56  if (x < xmin) xmin = x;
57  if (y < ymin) ymin = y;
58  if (z < zmin) zmin = z;
59  if (x > xmax) xmax = x;
60  if (y > ymax) ymax = y;
61  if (z > zmax) zmax = z;
62  }
63  delta_x = xmax - xmin + 1;
64  delta_y = ymax - ymin + 1;
65  delta_z = zmax - zmin + 1;
66 
67  // std::cout << "delta_x = " << delta_x
68  // << " delta_y = " << delta_y
69  // << " delta_z = " << delta_z << std::endl;
70  }
71 
72  // Compute dimension to split
74  split_dim = 0;
75  if (delta_y >= delta_x && delta_y >= delta_z) split_dim = 1;
76  if (delta_z >= delta_x && delta_z >= delta_y) split_dim = 2;
77 
78  //std::cout << "splitting dimension = " << split_dim << std::endl;
79  }
80 
81  // Split box into two pieces with roughly equal numbers of points
82  void split() {
83  // Sort points based on splitting dimension
85  std::sort(coords.begin(), coords.end(), cmp);
86 
87  // Divide coords into two bins of roughly equal size, keeping
88  // coords with equal values for split dimension together
89  size_type n = coords.size();
90  size_type s = n / 2;
91 
92  while (s < n-1 && coords[s](split_dim) == coords[s+1](split_dim)) ++s;
93  //std::cout << "n = " << n << " s = " << s << std::endl;
94 
95  if (s > 0)
96  left = Teuchos::rcp(new Box(coords.view(0, s)));
97  if (s < n)
98  right = Teuchos::rcp(new Box(coords.view(s, n-s)));
99 
100  // Clear my coordinate array since we aren't a leaf
101  //Teuchos::Array<TupleType>().swap(coords);
102  //coords.resize(0);
103  }
104 
105  };
106 
108  RCB(const coord_type& max_length_,
109  const size_type& max_parts_,
110  const Teuchos::ArrayView<TupleType>& coords_) :
111  max_length(max_length_),
112  max_parts(max_parts_),
113  coords(coords_.begin(), coords_.end()) {
114  partition();
115  }
116 
118  ~RCB() {}
119 
121  size_type get_num_parts() const { return num_parts; }
122 
125 
128  return parts; }
129 
130  // Create list of part IDs for each tuple
133  for (size_type part=0; part<num_parts; ++part) {
134  Teuchos::RCP<Box> box = (*parts)[part];
135  size_type n = box->coords.size();
136  for (size_type i=0; i<n; ++i)
137  part_ids[ box->coords[i].ID() ] = part;
138  }
139  return part_ids;
140  }
141 
142  private:
143 
149 
151  void partition() {
152 
153  // Create root bounding box
154  root = Teuchos::rcp(new Box(coords()));
155  num_parts = 1;
157 
158  // Create array of boxes that are too big
159  Teuchos::Array< Teuchos::RCP<Box> > boxes_to_split;
160  if (root->delta_x > max_length ||
161  root->delta_y > max_length ||
162  root->delta_z > max_length)
163  boxes_to_split.push_back(root);
164  else
165  parts->push_back(root);
166 
167  // Split each box until all boxes are less than tolerance
168  while(boxes_to_split.size() > 0 && num_parts < max_parts) {
169  Teuchos::RCP<Box> box = boxes_to_split.back();
170  boxes_to_split.pop_back();
171  box->split();
172  ++num_parts;
173  if (box->left != Teuchos::null) {
174  if (box->left->delta_x > max_length ||
175  box->left->delta_y > max_length ||
176  box->left->delta_z > max_length)
177  boxes_to_split.push_back(box->left);
178  else
179  parts->push_back(box->left);
180  }
181  if (box->right != Teuchos::null) {
182  if (box->right->delta_x > max_length ||
183  box->right->delta_y > max_length ||
184  box->right->delta_z > max_length)
185  boxes_to_split.push_back(box->right);
186  else
187  parts->push_back(box->right);
188  }
189  }
190 
191  TEUCHOS_ASSERT(parts->size() == num_parts);
192  }
193 
194  };
195 
196  template <typename ordinal_type, typename scalar_type>
197  struct CijkData {
203 
205  if (d == 0) return i;
206  if (d == 1) return j;
207  if (d == 2) return k;
208  return -1;
209  }
210 
211  ordinal_type ID() const { return gid; }
212  };
213 
218  };
219 
220  template <typename ordinal_type, typename scalar_type>
224  CijkSymmetryType symmetry_type) {
226  typedef typename Cijk_type::k_iterator k_iterator;
227  typedef typename Cijk_type::kj_iterator kj_iterator;
228  typedef typename Cijk_type::kji_iterator kji_iterator;
229 
230  ordinal_type num_cijk = Cijk.num_entries();
232  num_cijk);
233  ordinal_type idx = 0;
234  k_iterator k_begin = Cijk.k_begin();
235  k_iterator k_end = Cijk.k_end();
236  for (k_iterator k_it=k_begin; k_it!=k_end; ++k_it) {
237  ordinal_type k = index(k_it);
238  kj_iterator j_begin = Cijk.j_begin(k_it);
239  kj_iterator j_end = Cijk.j_end(k_it);
240  for (kj_iterator j_it = j_begin; j_it != j_end; ++j_it) {
241  ordinal_type j = index(j_it);
242  kji_iterator i_begin = Cijk.i_begin(j_it);
243  kji_iterator i_end = Cijk.i_end(j_it);
244  for (kji_iterator i_it = i_begin; i_it != i_end; ++i_it) {
245  ordinal_type i = index(i_it);
246  if (symmetry_type == CIJK_NO_SYMMETRY) {
247  coordinate_list[idx].i = i;
248  coordinate_list[idx].j = j;
249  coordinate_list[idx].k = k;
250  coordinate_list[idx].c = value(i_it);
251  coordinate_list[idx].gid = idx;
252  ++idx;
253  }
254  else if (symmetry_type == CIJK_TWO_WAY_SYMMETRY && j >= k) {
255  coordinate_list[idx].i = i;
256  coordinate_list[idx].j = j;
257  coordinate_list[idx].k = k;
258  if (j == k)
259  coordinate_list[idx].c = 0.5*value(i_it);
260  else
261  coordinate_list[idx].c = value(i_it);
262  coordinate_list[idx].gid = idx;
263  ++idx;
264  }
265  else if (symmetry_type == CIJK_SIX_WAY_SYMMETRY && i >= j && j >= k) {
266  coordinate_list[idx].i = i;
267  coordinate_list[idx].j = j;
268  coordinate_list[idx].k = k;
269  if (i == j && j == k)
270  coordinate_list[idx].c = (1.0/6.0)*value(i_it);
271  else
272  coordinate_list[idx].c = value(i_it);
273  coordinate_list[idx].gid = idx;
274  ++idx;
275  }
276  }
277  }
278  }
279  coordinate_list.resize(idx);
280 
281  return coordinate_list;
282  }
283 
284 } // namespace Stokhos
285 
286 #endif // STOKHOS_SPARSE_3_TENSOR_PARTITION_HPP
k_iterator k_begin() const
Iterator pointing to first k entry.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
Teuchos::Array< TupleType > coords
ArrayView< T > view(size_type offset, size_type size)
kj_iterator j_begin(const k_iterator &k) const
Iterator pointing to first j entry for given k.
size_type get_num_parts() const
Get number of parts.
Teuchos::ArrayRCP< id_type > get_part_IDs() const
kj_iterator j_end(const k_iterator &k) const
Iterator pointing to last j entry for given k.
void resize(const size_type n, const T &val=T())
Teuchos::RCP< Box > get_partition_root() const
Get root of partition.
Teuchos::RCP< Teuchos::Array< Teuchos::RCP< Box > > > parts
TupleType::value_type coord_type
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Box(const Teuchos::ArrayView< TupleType > &c)
bool operator()(const TupleType &a, const TupleType &b) const
Teuchos::Array< TupleType > coords
iterator end()
k_iterator k_end() const
Iterator pointing to last k entry.
Stokhos::Sparse3Tensor< int, double > Cijk_type
Teuchos::RCP< Teuchos::Array< Teuchos::RCP< Box > > > get_parts() const
Get parts array.
reference back()
void push_back(const value_type &x)
size_type size() const
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::ArrayView< TupleType >::size_type size_type
Teuchos::ArrayRCP< CijkData< ordinal_type, scalar_type > > build_cijk_coordinate_list(const Sparse3Tensor< ordinal_type, scalar_type > &Cijk, CijkSymmetryType symmetry_type)
kji_iterator i_begin(const kj_iterator &j) const
Iterator pointing to first i entry for given j and k.
kji_iterator i_end(const kj_iterator &j) const
Iterator pointing to last i entry for given j and k.
RCB(const coord_type &max_length_, const size_type &max_parts_, const Teuchos::ArrayView< TupleType > &coords_)
Constructor.
iterator begin()
int n
ordinal_type operator()(ordinal_type d) const
ordinal_type num_entries() const
Return number of non-zero entries.