Zoltan2
Loading...
Searching...
No Matches
vecWithCopies.cpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
16
17#include "Teuchos_CommHelpers.hpp"
18#include "Teuchos_DefaultComm.hpp"
19#include "Teuchos_RCP.hpp"
20#include "Teuchos_Array.hpp"
21#include "Tpetra_Core.hpp"
22#include "Tpetra_Map.hpp"
23#include "Tpetra_Vector.hpp"
24
25#include <string>
26#include <sstream>
27#include <iostream>
28
30
31int main(int narg, char **arg)
32{
33 typedef Tpetra::Map<> map_t;
34 typedef map_t::local_ordinal_type lno_t;
35 typedef map_t::global_ordinal_type gno_t;
36 typedef int scalar_t;
37
38 Tpetra::ScopeGuard tscope(&narg, &arg);
39 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
40 int me = comm->getRank();
41
42 // Create a map with duplicated entries (mapWithCopies)
43 // Each rank has 15 IDs, the last five of which overlap with the next rank.
44
45
46 lno_t numLocalCoords = 15;
47 lno_t offset = me * 10;
48
49 Teuchos::Array<gno_t> gids(numLocalCoords);
50 for (lno_t i = 0 ; i < numLocalCoords; i++)
51 gids[i] = static_cast<gno_t> (offset + i);
52
53 Tpetra::global_size_t numGlobalCoords =
54 Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid();
55 Teuchos::RCP<const map_t> mapWithCopies =
56 rcp(new map_t(numGlobalCoords, gids(), 0, comm));
57
58 // Create a new map with IDs uniquely assigned to ranks (oneToOneMap)
59 Teuchos::RCP<const map_t> oneToOneMap =
60 Tpetra::createOneToOne<lno_t, gno_t>(mapWithCopies);
61
62 // Create vectors with each map
63 typedef Tpetra::Vector<scalar_t, lno_t, gno_t> vector_t;
64
65 vector_t vecWithCopies(mapWithCopies);
66 vector_t oneToOneVec(oneToOneMap);
67
68 // Set values in oneToOneVec: each entry == rank
69 for (lno_t i = 0; i < lno_t(oneToOneMap->getLocalNumElements()); i++)
70 oneToOneVec.replaceLocalValue(i, me);
71
72 // Now import oneToOneVec's values back to vecWithCopies
73 Teuchos::RCP<const Tpetra::Import<lno_t, gno_t> > importer =
74 Tpetra::createImport<lno_t, gno_t>(oneToOneMap, mapWithCopies);
75 vecWithCopies.doImport(oneToOneVec, *importer, Tpetra::REPLACE);
76
77 // Print the entries of each vector
78 std::cout << me << " ONE TO ONE VEC ("
79 << oneToOneMap->getGlobalNumElements() << "): ";
80 lno_t nlocal = lno_t(oneToOneMap->getLocalNumElements());
81 for (lno_t i = 0; i < nlocal; i++)
82 std::cout << "[" << oneToOneMap->getGlobalElement(i) << " "
83 << oneToOneVec.getData()[i] << "] ";
84 std::cout << std::endl;
85
86 // Should see copied vector values when print VEC WITH COPIES
87 std::cout << me << " VEC WITH COPIES ("
88 << mapWithCopies->getGlobalNumElements() << "): ";
89 nlocal = lno_t(mapWithCopies->getLocalNumElements());
90 for (lno_t i = 0; i < nlocal; i++)
91 std::cout << "[" << mapWithCopies->getGlobalElement(i) << " "
92 << vecWithCopies.getData()[i] << "] ";
93 std::cout << std::endl;
94
95 return 0;
96}
int main()
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t
Tpetra::Map map_t