51 #include <Teuchos_DefaultComm.hpp>
52 #include <Teuchos_RCP.hpp>
53 #include <Teuchos_CommHelpers.hpp>
64 bool strideOne =
false;
66 if (valueStrides == NULL) strideOne =
true;
79 for (
int i=0; !fail && i < len; i++)
80 if (!fail && idList[i] != ids[i])
83 for (
int v=0; !fail && v < mvdim; v++){
85 int correctStride = (strideOne ? 1 : valueStrides[v]);
90 if (!fail && stride != correctStride)
93 for (
int i=0; !fail && i < len; i++){
100 for (
int w=0; !fail && w < wdim; w++){
106 if (!fail && stride != weightStrides[w])
109 for (
int i=0; !fail && i < len; i++){
110 if (wgts[stride*i] != weights[w][weightStrides[w]*i])
119 int main(
int narg,
char *arg[])
121 Tpetra::ScopeGuard tscope(&narg, &arg);
122 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
124 int rank = comm->getRank();
125 int nprocs = comm->getSize();
133 zgno_t myFirstId = rank * numLocalIds;
137 int *weightStrides =
new int [wdim];
143 int *valueStrides =
new int [mvdim];
146 for (
zlno_t i=0; i < numLocalIds; i++){
147 myIds[i] = myFirstId+i;
149 for (
int w=0; w < wdim; w++)
150 weights[w*numLocalIds + i] = w + 1 + nprocs - rank;
152 v_values[i] = numLocalIds-i;
154 for (
int v=0; v < mvdim; v++)
155 mv_values[i*mvdim + v] = (v+1) * (nprocs-rank) / (i+1);
158 for (
int w=0; w < wdim; w++){
159 weightStrides[w] = 1;
160 weightPtrs[w] = weights + numLocalIds*w;
163 for (
int v=0; v < mvdim; v++){
164 valueStrides[v] = mvdim;
165 valuePtrs[v] = mv_values + v;
173 std::vector<const zscalar_t *> weightValues;
174 std::vector<int> strides;
180 catch (std::exception &e){
187 myIds, 1, valuePtrs, NULL, 0, NULL, NULL);
197 std::vector<const zscalar_t *> weightValues;
198 std::vector<int> strides;
200 weightValues.push_back(weightPtrs[0]);
201 weightValues.push_back(weightPtrs[1]);
202 strides.push_back(1);
203 strides.push_back(1);
207 v_values, 1,
true, weightPtrs[0], 1);
209 catch (std::exception &e){
216 myIds, 1, valuePtrs, NULL, 1, weightPtrs, weightStrides);
226 std::vector<const zscalar_t *> weightValues, values;
227 std::vector<int> wstrides, vstrides;
229 for (
int dim=0; dim < mvdim; dim++){
230 values.push_back(valuePtrs[dim]);
231 vstrides.push_back(mvdim);
237 numLocalIds, myIds, values, vstrides, weightValues, wstrides);
239 catch (std::exception &e){
246 myIds, mvdim, valuePtrs, valueStrides, 0, NULL, NULL);
256 std::vector<const zscalar_t *> weightValues, values;
257 std::vector<int> wstrides, vstrides;
259 for (
int dim=0; dim < wdim; dim++){
260 weightValues.push_back(weightPtrs[dim]);
261 wstrides.push_back(1);
264 for (
int dim=0; dim < mvdim; dim++){
265 values.push_back(valuePtrs[dim]);
266 vstrides.push_back(mvdim);
271 numLocalIds, myIds, values, vstrides, weightValues, wstrides);
274 catch (std::exception &e){
281 myIds, mvdim, valuePtrs, valueStrides,
282 wdim, weightPtrs, weightStrides);
290 std::cout <<
"PASS" << std::endl;
294 delete [] weightStrides;
295 delete [] weightPtrs;
298 delete [] valueStrides;
size_t getLocalNumIDs() const
Returns the number of objects on this process.
int main(int narg, char *arg[])
common code used by tests
int getNumEntriesPerID() const
Return the number of vectors.
list idList
Match up parameters to validators.
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
static const std::string fail
#define TEST_FAIL_AND_RETURN_VALUE(comm, ok, s, rc)
void getIDsView(const gno_t *&ids) const
Provide a pointer to this process' identifiers.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater...
void getEntriesView(const scalar_t *&entries, int &stride, int idx=0) const
Defines the BasicVectorAdapter class.