Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_SolverMap_CrsMatrix_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef TPETRA_SOLVERMAP_CRSMATRIX_DEF_HPP
11 #define TPETRA_SOLVERMAP_CRSMATRIX_DEF_HPP
12 
22 
24 
25 #include <vector>
26 #include <filesystem>
27 
28 namespace Tpetra {
29 
30 template <class Scalar,
31  class LocalOrdinal,
32  class GlobalOrdinal,
33  class Node>
35  : StructuralSameTypeTransform<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >()
36  , newColMap_(Teuchos::null)
37  , newGraph_(Teuchos::null) {
38  // Nothing to do
39 }
40 
41 template <class Scalar,
42  class LocalOrdinal,
43  class GlobalOrdinal,
44  class Node>
46  // Nothing to do
47 }
48 
49 template <class Scalar,
50  class LocalOrdinal,
51  class GlobalOrdinal,
52  class Node>
53 typename SolverMap_CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::NewType
56 
57  assert(!origMatrix->isGloballyIndexed());
58 
59  this->origObj_ = origMatrix;
60 
61  // *******************************************************************
62  // Step 1/7: Check if domain map and col map are different
63  // *******************************************************************
64  Teuchos::RCP<map_t const> origDomainMap = origMatrix->getDomainMap();
65  Teuchos::RCP<map_t const> origColMap = origMatrix->getColMap();
66 
67  typename map_t::local_map_type localOrigDomainMap(origDomainMap->getLocalMap());
68  typename map_t::local_map_type localOrigColMap(origColMap->getLocalMap());
69 
70  if (origDomainMap->isLocallyFitted(*origColMap)) {
71  this->newObj_ = this->origObj_;
72  } else {
75 
76  // *****************************************************************
77  // Step 2/7: Fill newColMap_globalColIndices with the global indices
78  // of all entries in origDomainMap
79  // *****************************************************************
80 
81  // Initialize the value of 'newColMap_localSize'
82  size_t const origDomainMap_localSize = origDomainMap->getLocalNumElements();
83  size_t newColMap_localSize(origDomainMap_localSize);
84 
85  // Increase the value of 'newColMap_localSize' as necessary
86  size_t newColMap_extraSize(0);
87  size_t const origColMap_localSize = origColMap->getLocalNumElements();
88  {
89  Kokkos::parallel_reduce(
90  "Tpetra::SolverMap_CrsMatrix::construct::newColMap_extraSize", Kokkos::RangePolicy<typename Node::device_type::execution_space, size_t>(0, origColMap_localSize), KOKKOS_LAMBDA(size_t const i, size_t& sizeToUpdate)->void {
91  GlobalOrdinal const globalColIndex(localOrigColMap.getGlobalElement(i));
92  if (localOrigDomainMap.getLocalElement(globalColIndex) == ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid()) {
93  sizeToUpdate += 1;
94  }
95  },
96  newColMap_extraSize);
97  }
98  newColMap_localSize += newColMap_extraSize;
99 
100  // Instantiate newColMap_globalColIndices with the correct size 'newColMap_localSize'
101  Kokkos::View<GlobalOrdinal*, typename Node::device_type> newColMap_globalColIndices("", newColMap_localSize);
102 
103  // Fill newColMap_globalColIndices with the global indices of all entries in origDomainMap
104  {
105  Kokkos::parallel_for(
106  "Tpetra::SolverMap_CrsMatrix::construct::copyDomainMapToNewColMap", Kokkos::RangePolicy<typename Node::device_type::execution_space, size_t>(0, origDomainMap_localSize), KOKKOS_LAMBDA(size_t const i)->void {
107  newColMap_globalColIndices(i) = localOrigDomainMap.getGlobalElement(i);
108  });
109  }
110 
111  // *****************************************************************
112  // Step 3/7: Append newColMap_globalColIndices with those origColMap
113  // entries that are not in newColMap_globalColIndices yet
114  // *****************************************************************
115  {
116  Kokkos::parallel_scan(
117  "Tpetra::SolverMap_CrsMatrix::construct::appendNewColMap", Kokkos::RangePolicy<typename Node::device_type::execution_space, size_t>(0, origColMap_localSize), KOKKOS_LAMBDA(size_t const i, size_t& jToUpdate, bool const final)->void {
118  GlobalOrdinal const globalColIndex(localOrigColMap.getGlobalElement(i));
119  if (localOrigDomainMap.getLocalElement(globalColIndex) == ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid()) {
120  if (final) {
121  newColMap_globalColIndices(origDomainMap_localSize + jToUpdate) = globalColIndex;
122  }
123  jToUpdate += 1;
124  }
125  });
126  }
127 
128  // *****************************************************************
129  // Step 4/7: Create a new column map using newColMap_globalColIndices
130  // *****************************************************************
131  Teuchos::RCP<map_t const> origRowMap = origMatrix->getRowMap();
132  Teuchos::RCP<Teuchos::Comm<int> const> Comm = origRowMap->getComm();
133  size_t const newColMap_localNumCols = newColMap_globalColIndices.size();
134  size_t newColMap_globalNumCols(0);
135  Teuchos::reduceAll(*Comm, Teuchos::REDUCE_SUM, 1, &newColMap_localNumCols, &newColMap_globalNumCols);
136 
137  newColMap_ = Teuchos::rcp<map_t>(new map_t(newColMap_globalNumCols, newColMap_globalColIndices, origDomainMap->getIndexBase(), Comm));
138 
139  // *****************************************************************
140  // Step 5/7: Create new graph
141  // *****************************************************************
142  cg_t const* origGraph = dynamic_cast<cg_t const*>(origMatrix->getGraph().get());
143  newGraph_ = Teuchos::rcp<cg_t>(new cg_t(*origGraph));
144  newGraph_->resumeFill();
145  newGraph_->reindexColumns(newColMap_);
146  newGraph_->fillComplete();
147 
148  // *****************************************************************
149  // Step 6/7: Create new CRS matrix
150  // *****************************************************************
151  Teuchos::RCP<cm_t> newMatrix = Teuchos::rcp<cm_t>(new cm_t(*origMatrix, newGraph_));
152 
153  // *****************************************************************
154  // Step 7/7: Update newObj_
155  // *****************************************************************
156  this->newObj_ = newMatrix;
157  }
158 
159  return this->newObj_;
160 }
161 
162 //
163 // Explicit instantiation macro
164 //
165 // Must be expanded from within the Tpetra namespace!
166 //
167 
168 #define TPETRA_SOLVERMAPCRSMATRIX_INSTANT(SCALAR, LO, GO, NODE) \
169  template class SolverMap_CrsMatrix<SCALAR, LO, GO, NODE>;
170 
171 } // namespace Tpetra
172 
173 #endif // TPETRA_SOLVERMAP_CRSMATRIX_DEF_HPP
Sparse matrix that presents a row-oriented interface that lets users read or modify entries...
Declaration of the Tpetra::SolverMap_CrsMatrix class.
NewType operator()(OriginalType const &origMatrix)
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
A parallel distribution of indices over processes.