Zoltan2
Toggle main menu visibility
Loading...
Searching...
No Matches
Zoltan2_StridedData.hpp
Go to the documentation of this file.
1
// @HEADER
2
//
3
// ***********************************************************************
4
//
5
// Zoltan2: A package of combinatorial algorithms for scientific computing
6
// Copyright 2012 Sandia Corporation
7
//
8
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9
// the U.S. Government retains certain rights in this software.
10
//
11
// Redistribution and use in source and binary forms, with or without
12
// modification, are permitted provided that the following conditions are
13
// met:
14
//
15
// 1. Redistributions of source code must retain the above copyright
16
// notice, this list of conditions and the following disclaimer.
17
//
18
// 2. Redistributions in binary form must reproduce the above copyright
19
// notice, this list of conditions and the following disclaimer in the
20
// documentation and/or other materials provided with the distribution.
21
//
22
// 3. Neither the name of the Corporation nor the names of the
23
// contributors may be used to endorse or promote products derived from
24
// this software without specific prior written permission.
25
//
26
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37
//
38
// Questions? Contact Karen Devine (kddevin@sandia.gov)
39
// Erik Boman (egboman@sandia.gov)
40
// Siva Rajamanickam (srajama@sandia.gov)
41
//
42
// ***********************************************************************
43
//
44
// @HEADER
45
#ifndef _ZOLTAN2_STRIDEDDATA_HPP_
46
#define _ZOLTAN2_STRIDEDDATA_HPP_
47
48
#include <
Zoltan2_Standards.hpp
>
49
#include <
Zoltan2_Environment.hpp
>
50
#include <typeinfo>
51
55
56
namespace
Zoltan2
{
57
74
75
template
<
typename
lno_t,
typename
scalar_t>
76
class
StridedData
{
77
private
:
78
ArrayRCP<const scalar_t> vec_;
79
int
stride_;
80
81
public
:
82
89
StridedData
(ArrayRCP<const scalar_t> x,
int
stride) :
90
vec_(x), stride_(stride)
91
{ }
92
95
StridedData
(): vec_(), stride_(0)
96
{ }
97
103
lno_t
size
()
const
{
return
vec_.size(); }
104
111
scalar_t
operator[]
(
lno_t
idx)
const
{
return
vec_[idx*stride_]; }
112
121
122
template
<
typename
T>
void
getInputArray
(ArrayRCP<const T> &array)
const
;
123
129
void
getStridedList
(ArrayRCP<const scalar_t> &vec,
int
&stride)
const
130
{
131
vec = vec_;
132
stride = stride_;
133
}
134
142
void
getStridedList
(
size_t
&len,
const
scalar_t
*&vec,
int
&stride)
const
143
{
144
len = vec_.size();
145
if
(len != 0) vec = vec_.getRawPtr();
146
else
vec = NULL;
147
stride = stride_;
148
}
149
152
StridedData
&
operator=
(
const
StridedData
&sInput)
153
{
154
if
(
this
!= &sInput)
155
sInput.
getStridedList
(vec_, stride_);
156
157
return
*
this
;
158
}
159
};
160
161
// Helper function needed for T=scalar_t specialization
162
// Separate function needed // because, with T != scalar_t,
163
// "array = vec_" // would not compile;
164
// ArrayRCP does not overload "=" operator for different types.
165
// Separate helper function needed (outside StridedData class)
166
// because cannot specialize member function of templated class.
167
template
<
typename
scalar_t,
typename
T>
168
static
void
getInputArrayHelper
(
169
ArrayRCP<const T> &target,
170
const
ArrayRCP<const scalar_t> &src)
171
{
172
// Create a copy of desired type T
173
// From logic in getInputArray, we know stride == 1.
174
size_t
n = src.size();
175
T *tmp =
new
T [n];
176
177
if
(!tmp){
178
std::cerr <<
"Error: "
<< __FILE__ <<
", "
<< __LINE__<< std::endl;
179
std::cerr << n <<
" objects"
<< std::endl;
180
throw
std::bad_alloc();
181
}
182
183
for
(
size_t
i=0; i < n; i++){
184
tmp[i] =
static_cast<
T
>
(src[i]);
185
}
186
target = arcp(tmp, 0, n);
187
}
188
189
// Specialization with T == scalar_t: just copy ArrayRCP
190
template
<
typename
scalar_t>
191
static
void
getInputArrayHelper
(
192
ArrayRCP<const scalar_t> &target,
193
const
ArrayRCP<const scalar_t> &src)
194
{
195
target = src;
196
}
197
198
// Member function for getting unstrided view/copy of StridedData.
199
template
<
typename
lno_t,
typename
scalar_t>
200
template
<
typename
T>
201
void
StridedData<lno_t, scalar_t>::getInputArray
(
202
ArrayRCP<const T> &array)
const
203
{
204
if
(vec_.size() < 1){
205
array = ArrayRCP<const T>();
206
}
207
else
if
(stride_ > 1) {
208
// Create an unstrided copy
209
size_t
n = vec_.size() / stride_;
210
T *tmp =
new
T [n];
211
212
if
(!tmp){
213
std::cerr <<
"Error: "
<< __FILE__ <<
", "
<< __LINE__<< std::endl;
214
std::cerr << n <<
" objects"
<< std::endl;
215
throw
std::bad_alloc();
216
}
217
218
for
(
size_t
i=0,j=0; i < n; i++,j+=stride_){
219
tmp[i] =
static_cast<
T
>
(vec_[j]);
220
}
221
array = arcp(tmp, 0, n);
222
}
223
else
{
// stride == 1
224
Zoltan2::getInputArrayHelper<scalar_t, T>
(array, vec_);
225
}
226
return
;
227
}
228
229
}
// namespace Zoltan2
230
231
#endif
Zoltan2_Environment.hpp
Defines the Environment class.
Zoltan2_Standards.hpp
Gathering definitions used in software development.
StridedData::StridedData
StridedData(ArrayRCP< const scalar_t > x, int stride)
Constructor.
Definition
Zoltan2_StridedData.hpp:89
Zoltan2::StridedData::getStridedList
void getStridedList(size_t &len, const scalar_t *&vec, int &stride) const
Get the raw input information.
Definition
Zoltan2_StridedData.hpp:142
Zoltan2::StridedData::StridedData
StridedData(ArrayRCP< const scalar_t > x, int stride)
Constructor.
Definition
Zoltan2_StridedData.hpp:89
Zoltan2::StridedData::operator=
StridedData & operator=(const StridedData &sInput)
Assignment operator.
Definition
Zoltan2_StridedData.hpp:152
Zoltan2::StridedData::getInputArray
void getInputArray(ArrayRCP< const T > &array) const
Create a contiguous array of the required type, perhaps for a TPL.
Definition
Zoltan2_StridedData.hpp:201
Zoltan2::StridedData::operator[]
scalar_t operator[](lno_t idx) const
Access an element of the input array.
Definition
Zoltan2_StridedData.hpp:111
Zoltan2::StridedData::getStridedList
void getStridedList(ArrayRCP< const scalar_t > &vec, int &stride) const
Get a reference counted pointer to the input.
Definition
Zoltan2_StridedData.hpp:129
Zoltan2::StridedData::size
lno_t size() const
Return the length of the strided array.
Definition
Zoltan2_StridedData.hpp:103
Zoltan2::StridedData::StridedData
StridedData()
Default constructor. A zero-length strided array.
Definition
Zoltan2_StridedData.hpp:95
Zoltan2
Created by mbenlioglu on Aug 31, 2020.
Definition
Zoltan2_AlgHybrid2GL.hpp:38
Zoltan2::getInputArrayHelper
static void getInputArrayHelper(ArrayRCP< const T > &target, const ArrayRCP< const scalar_t > &src)
Definition
Zoltan2_StridedData.hpp:168
Zoltan2::lno_t
core
src
util
Zoltan2_StridedData.hpp
Generated by
1.17.0