File: TpetraUtils_crsUtils.cpp

package info (click to toggle)
trilinos 12.18.1-2
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 825,604 kB
  • sloc: cpp: 3,352,065; ansic: 432,926; fortran: 169,495; xml: 61,903; python: 40,836; sh: 32,697; makefile: 26,612; javascript: 8,535; perl: 7,136; f90: 6,372; csh: 4,183; objc: 2,620; lex: 1,469; lisp: 810; yacc: 497; awk: 364; ml: 281; php: 145
file content (341 lines) | stat: -rw-r--r-- 11,875 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
/*
// @HEADER
// ***********************************************************************
//
//          Tpetra: Templated Linear Algebra Services Package
//                 Copyright (2008) Sandia Corporation
//
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
//
// ************************************************************************
// @HEADER
*/

#include "Tpetra_TestingUtilities.hpp"
#include "Tpetra_Core.hpp"
#include "Tpetra_Map.hpp"
#include "Tpetra_ConfigDefs.hpp"
#include "Tpetra_Details_crsUtils.hpp"

#include "Kokkos_Core.hpp"
#include "Teuchos_UnitTestHarness.hpp"
#include "Teuchos_Array.hpp"
#include "Teuchos_as.hpp"

#include <algorithm>
#include <iterator>

using Teuchos::CommandLineProcessor;
using Tpetra::Details::padCrsArrays;
using Tpetra::Details::insertCrsIndices;
using Tpetra::Details::findCrsIndices;
using Tpetra::Details::impl::uninitialized_view;
using std::vector;

namespace {

  TEUCHOS_STATIC_SETUP()
  {
    Teuchos::CommandLineProcessor &clp = Teuchos::UnitTestRepository::getCLP();
    clp.addOutputSetupOptions(true);
  }

  //
  // UNIT TESTS
  //
TEUCHOS_UNIT_TEST(CrsGraph, ResizeRowPointersAndIndices_1)
{
  using device_type = typename Tpetra::Map<>::device_type;
  using execution_space = typename device_type::execution_space;
  using ordinal_type = size_t;
  using view_type = Kokkos::View<ordinal_type*, device_type>;
  using size_type = typename view_type::size_type;

  ordinal_type num_row = 4;
  ordinal_type num_indices_per_row = 5;
  ordinal_type num_indices = num_indices_per_row * num_row;
  auto row_ptrs_beg = uninitialized_view<view_type>("beg", num_row+1);
  // this assumes UVM
  for (ordinal_type i=0; i<num_row+1; i++) row_ptrs_beg(i) = num_indices_per_row*i;

  auto row_ptrs_end = uninitialized_view<view_type>("end", num_row);
  for (ordinal_type i=0; i<num_row; i++) row_ptrs_end(i) = row_ptrs_beg(i+1);

  auto indices = uninitialized_view<view_type>("indices", num_indices);
  for (ordinal_type i=0; i<num_row; i++) {
    auto start = row_ptrs_beg(i);
    auto end = row_ptrs_beg(i+1);
    for (ordinal_type j=start, k=0; j<end; j++, k++) {
      indices(j) = (i + 1) * num_indices_per_row + k;
    }
  }

  auto import_lids = uninitialized_view<view_type>("import lids", num_row);
  auto num_packets_per_lid = uninitialized_view<view_type>("num packets", num_row);
  for (ordinal_type i=0; i<num_row; i++) {
   import_lids(i) = i;
   num_packets_per_lid(i) = i;
  }
  ordinal_type num_extra =
    num_row*(num_packets_per_lid(0) + num_packets_per_lid(num_row-1))/2;

  Kokkos::UnorderedMap<ordinal_type,ordinal_type,device_type> padding(import_lids.size());
  execution_space().fence();
  for (size_type i=0; i<import_lids.size(); i++){
    padding.insert(import_lids(i), num_packets_per_lid(i));
  }
  execution_space().fence();
  TEST_ASSERT(!padding.failed_insert());

  padCrsArrays(row_ptrs_beg, row_ptrs_end, indices, padding);
  TEST_ASSERT(indices.size() == static_cast<size_type>(num_indices + num_extra));

  {
    // make sure row_ptrs_beg is where it should be
    bool offsets_ok = row_ptrs_beg(0) == 0;
    for (ordinal_type i=1; i<num_row; i++) {
      auto expected = row_ptrs_beg(i-1) + num_indices_per_row + num_packets_per_lid(i-1);
      if (row_ptrs_beg(i) != expected) offsets_ok = false;
    }
    if (offsets_ok) offsets_ok = row_ptrs_beg(num_row) == indices.size();
    TEST_ASSERT(offsets_ok);
  }

  {
    // make sure indices were shifted correctly
    bool indices_ok = true;
    for (ordinal_type i=0; i<num_row; i++) {
      auto start = row_ptrs_beg(i);
      auto end = row_ptrs_end(i);
      for (ordinal_type j=start, k=0; j<end; j++, k++) {
        auto expected = (i + 1) * num_indices_per_row + k;
        if (expected != indices(j)) indices_ok = false;
      }
    }
    TEST_ASSERT(indices_ok);
  }
}

TEUCHOS_UNIT_TEST(CrsGraph, ResizeRowPointersAndIndices_2)
{
  typedef typename Tpetra::Map<>::device_type device_type;
  using execution_space = typename device_type::execution_space;
  using ordinal_type = size_t;
  using view_type = Kokkos::View<ordinal_type*, device_type>;
  using size_type = typename view_type::size_type;

  auto row_ptrs_beg = view_type("beg", 4);
  auto row_ptrs_end = view_type("end", 3);
  auto indices = view_type("indices", 9);

  // Row 1, 3 allocations, 3 used
  row_ptrs_beg(0) = 0; row_ptrs_end(0) = 3;
  indices(0) = 1; indices(1) = 2; indices(2) = 3;  // Row 1

  // Row 2, 3 allocations only 1 used
  row_ptrs_beg(1) = 3; row_ptrs_end(1) = 4;
  indices(3) = 4;

  // Row 3, 3 allocations, only 2 used
  row_ptrs_beg(2) = 6; row_ptrs_end(2) = 8;
  indices(6) = 7; indices(7) = 8;

  row_ptrs_beg(3) = 9;

  // Import 5 extra values in to Row 1 and 3 extra in to Row 3
  auto import_lids = view_type("import lids", 2);
  auto num_packets_per_lid = view_type("num packets", 2);

  // Import LIDs not ordered
  import_lids(0) = 2;
  num_packets_per_lid(0) = 3;

  import_lids(1) = 0;
  num_packets_per_lid(1) = 5;

  Kokkos::UnorderedMap<ordinal_type,ordinal_type,device_type> padding(import_lids.size());
  execution_space().fence();
  for (size_type i=0; i<import_lids.size(); i++){
    padding.insert(import_lids(i), num_packets_per_lid(i));
  }
  execution_space().fence();
  TEST_ASSERT(!padding.failed_insert());
  padCrsArrays(row_ptrs_beg, row_ptrs_end, indices, padding);

  // Check row offsets
  TEST_ASSERT(row_ptrs_beg(0) == 0);
  TEST_ASSERT(row_ptrs_beg(1) == 8);
  TEST_ASSERT(row_ptrs_beg(2) == 11);
  TEST_ASSERT(row_ptrs_beg(3) == 16);
  TEST_ASSERT(indices.size() == 16);

  TEST_ASSERT(row_ptrs_end(0) == 3);
  TEST_ASSERT(row_ptrs_end(1) == 9);
  TEST_ASSERT(row_ptrs_end(2) == 13);

  // Row 1
  TEST_ASSERT(indices(0) == 1);
  TEST_ASSERT(indices(1) == 2);
  TEST_ASSERT(indices(2) == 3);
  // 5 extra entries in Row 1

  // Row 2
  TEST_ASSERT(indices(8) == 4);
  // 2 unused indices in Row 2

  // Row 2
  TEST_ASSERT(indices(11) == 7);
  TEST_ASSERT(indices(12) == 8);

}

template <class V1, class V2>
bool
compare_array_values(V1 const& arr1, V2 const& arr2, size_t const n)
{
  bool l_success = true;
  for (size_t i = 0; i < n; i++)
  {
    if (arr1[i] != arr2[i])
    {
      l_success = false;
      //std::cout << "ARR1[" << i << "] = " << arr1[i] << ", expected " << arr2[i] << "\n";
    }
  }
  return l_success;
}

TEUCHOS_UNIT_TEST( TpetraUtils, insertIndices )
{
  {
    vector<int> row_ptrs{0, 10};
    vector<int> cur_indices{1, 3, 5, 7, -1, -1, -1, -1, -1, -1};
    size_t num_assigned = 4;
    vector<int> new_indices{0, 2, 4, 6, 8, 7, 5, 3, 1, 2, 6, 8, 4, 0};
    auto num_inserted = insertCrsIndices(0, row_ptrs, cur_indices, num_assigned, new_indices);
    TEST_EQUALITY(static_cast<int>(num_inserted), 5);
    vector<int> expected{1, 3, 5, 7, 0, 2, 4, 6, 8};
    TEST_ASSERT(compare_array_values(cur_indices, expected, expected.size()));
  }

  {
    vector<int> row_ptrs{0, 8};
    vector<int> cur_indices{3, 6, 9, 12, -1, -1, -1, -1};
    size_t num_assigned = 4;
    vector<int> new_indices{1, 0, 2};
    vector<int> expected{3, 6, 9, 12, 1, 0, 2, -1};
    auto num_inserted = insertCrsIndices(0, row_ptrs, cur_indices, num_assigned, new_indices);
    TEST_EQUALITY(static_cast<size_t>(num_inserted),
                  static_cast<size_t>(new_indices.size()));
    TEST_ASSERT(compare_array_values(cur_indices, expected, expected.size()));
  }

  {
    vector<int> row_ptrs{0, 7};
    vector<int> cur_indices{3, 6, 9, 12, -1, -1, -1};
    size_t num_assigned = 4;
    vector<int> new_indices{1, 0, 2, 5};
    auto num_inserted = insertCrsIndices(0, row_ptrs, cur_indices, num_assigned, new_indices);
    TEST_EQUALITY(static_cast<int>(num_inserted), -1);
  }

  {
    vector<int> row_ptrs{0, 9};
    vector<int> cur_indices{3, 6, 9, 12, -1, -1, -1, -1, -1};
    size_t num_assigned = 4;
    vector<int> new_indices{0, 4, 7, 10, 13};
    vector<int> expected{3, 6, 9, 12, 0, 4, 7, 10, 13};
    auto num_inserted = insertCrsIndices(0, row_ptrs, cur_indices, num_assigned, new_indices);
    TEST_EQUALITY(num_inserted, new_indices.size());
    TEST_ASSERT(compare_array_values(cur_indices, expected, expected.size()));
  }
}

TEUCHOS_UNIT_TEST( TpetraUtils, insertIndicesWithCallback )
{
  {
    vector<int> row_ptrs{0, 7};
    vector<int> cur_indices{3, 6, 9, 12, -1, -1, -1};
    size_t num_assigned = 4;
    vector<int> new_indices{1, 0, 2, 2, 1, 0, 1, 2, 1};
    vector<int> in_values{1, 1, 1, 1, 1, 1, 1, 1, 1};
    vector<int> expected{3, 6, 9, 12, 1, 0, 2};
    vector<int> values(cur_indices.size(), 0);
    vector<int> expected_values{0, 0, 0, 0, 4, 2, 3};
    auto num_inserted =
      insertCrsIndices(0, row_ptrs, cur_indices, num_assigned, new_indices,
        [&](const size_t k, const size_t start, const size_t offset){
          values[start+offset] += in_values[k];
        });
    for (size_t k=0; k<expected_values.size(); k++)
    {
      std::cout << "[" << k << "] = (" << expected_values[k] << ", " << values[k] << ")\n";
    }
    TEST_EQUALITY(num_inserted, 3);
    TEST_ASSERT(compare_array_values(cur_indices, expected, expected.size()));
    TEST_ASSERT(compare_array_values(values, expected_values, expected_values.size()));
  }
}

TEUCHOS_UNIT_TEST( TpetraUtils, findIndices )
{
  {
    vector<int> row_ptrs{0, 4};
    vector<int> cur_indices{3, 6, 9, 12};
    const size_t cur_num_entries = 4;
    vector<int> new_indices{3, 6, 9, 12, 12, 9, 3, 6, 0, 2};
    vector<int> in_values(new_indices.size(), 1);
    vector<int> values(cur_indices.size(), 0);
    vector<int> expected_values{2, 2, 2, 2};
    auto num_found =
      findCrsIndices(0, row_ptrs, cur_num_entries, cur_indices, new_indices,
        [&](const size_t k, const size_t start, const size_t offset){
          values[start+offset] += in_values[k];
        });
    TEST_EQUALITY(num_found, 8);
    TEST_ASSERT(compare_array_values(values, expected_values, expected_values.size()));
  }
}

} // namespace (anonymous)

int
main (int argc, char* argv[])
{
  // Initialize MPI (if enabled) before initializing Kokkos.  This
  // lets MPI control things like pinning processes to sockets.
  Tpetra::ScopeGuard tpetraScope (&argc, &argv);
  const int errCode =
    Teuchos::UnitTestRepository::runUnitTestsFromMain (argc, argv);
  return errCode;
}