42 #ifndef STOKHOS_SPARSE_3_TENSOR_PARTITION_HPP
43 #define STOKHOS_SPARSE_3_TENSOR_PARTITION_HPP
47 #include "Teuchos_ArrayRCP.hpp"
48 #include "Teuchos_Array.hpp"
52 template <
typename TupleType>
56 typedef typename Teuchos::ArrayView<TupleType>::size_type
size_type;
57 typedef typename TupleType::id_type
id_type;
62 bool operator() (
const TupleType& a,
const TupleType& b)
const {
74 Box(
const Teuchos::ArrayView<TupleType>& c) :
coords(c.begin(), c.end()) {
142 const Teuchos::ArrayView<TupleType>& coords_) :
145 coords(coords_.begin(), coords_.end()) {
159 Teuchos::RCP< Teuchos::Array< Teuchos::RCP<Box> > >
get_parts()
const {
164 Teuchos::ArrayRCP<id_type> part_ids(
coords.size());
166 Teuchos::RCP<Box> box = (*parts)[part];
169 part_ids[ box->coords[i].ID() ] = part;
180 Teuchos::RCP< Teuchos::Array< Teuchos::RCP<Box> > >
parts;
188 parts = Teuchos::rcp(
new Teuchos::Array< Teuchos::RCP<Box> >);
191 Teuchos::Array< Teuchos::RCP<Box> > boxes_to_split;
195 boxes_to_split.push_back(
root);
201 Teuchos::RCP<Box> box = boxes_to_split.back();
202 boxes_to_split.pop_back();
205 if (box->left != Teuchos::null) {
209 boxes_to_split.push_back(box->left);
211 parts->push_back(box->left);
213 if (box->right != Teuchos::null) {
217 boxes_to_split.push_back(box->right);
219 parts->push_back(box->right);
228 template <
typename ordinal_type,
typename scalar_type>
237 if (d == 0)
return i;
238 if (d == 1)
return j;
239 if (d == 2)
return k;
252 template <
typename ordinal_type,
typename scalar_type>
253 Teuchos::ArrayRCP< CijkData<ordinal_type,scalar_type> >
258 typedef typename Cijk_type::k_iterator k_iterator;
259 typedef typename Cijk_type::kj_iterator kj_iterator;
260 typedef typename Cijk_type::kji_iterator kji_iterator;
263 Teuchos::ArrayRCP< CijkData<ordinal_type,scalar_type> > coordinate_list(
266 k_iterator k_begin = Cijk.
k_begin();
267 k_iterator k_end = Cijk.
k_end();
268 for (k_iterator k_it=k_begin; k_it!=k_end; ++k_it) {
270 kj_iterator j_begin = Cijk.
j_begin(k_it);
271 kj_iterator j_end = Cijk.
j_end(k_it);
272 for (kj_iterator j_it = j_begin; j_it != j_end; ++j_it) {
274 kji_iterator i_begin = Cijk.
i_begin(j_it);
275 kji_iterator i_end = Cijk.
i_end(j_it);
276 for (kji_iterator i_it = i_begin; i_it != i_end; ++i_it) {
279 coordinate_list[idx].i = i;
280 coordinate_list[idx].j =
j;
281 coordinate_list[idx].k = k;
282 coordinate_list[idx].c = value(i_it);
283 coordinate_list[idx].gid = idx;
287 coordinate_list[idx].i = i;
288 coordinate_list[idx].j =
j;
289 coordinate_list[idx].k = k;
291 coordinate_list[idx].c = 0.5*value(i_it);
293 coordinate_list[idx].c = value(i_it);
294 coordinate_list[idx].gid = idx;
298 coordinate_list[idx].i = i;
299 coordinate_list[idx].j =
j;
300 coordinate_list[idx].k = k;
301 if (i ==
j &&
j == k)
302 coordinate_list[idx].c = (1.0/6.0)*value(i_it);
304 coordinate_list[idx].c = value(i_it);
305 coordinate_list[idx].gid = idx;
311 coordinate_list.resize(idx);
313 return coordinate_list;
318 #endif // STOKHOS_SPARSE_3_TENSOR_PARTITION_HPP