Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_EquilibrationInfo.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_DETAILS_EQUILIBRATIONINFO_HPP
11 #define TPETRA_DETAILS_EQUILIBRATIONINFO_HPP
12 
15 
16 #include "TpetraCore_config.h"
17 #if KOKKOS_VERSION >= 40799
18 #include "KokkosKernels_ArithTraits.hpp"
19 #else
20 #include "Kokkos_ArithTraits.hpp"
21 #endif
22 #include "Kokkos_Core.hpp"
23 
24 namespace Tpetra {
25 namespace Details {
26 
50 template <class ScalarType, class DeviceType>
52 #if KOKKOS_VERSION >= 40799
53  using val_type = typename KokkosKernels::ArithTraits<ScalarType>::val_type;
54 #else
55  using val_type = typename Kokkos::ArithTraits<ScalarType>::val_type;
56 #endif
57 #if KOKKOS_VERSION >= 40799
58  using mag_type = typename KokkosKernels::ArithTraits<val_type>::mag_type;
59 #else
60  using mag_type = typename Kokkos::ArithTraits<val_type>::mag_type;
61 #endif
62  using device_type = typename DeviceType::device_type;
63  using host_device_type = typename Kokkos::View<mag_type*, device_type>::host_mirror_type::device_type;
65 
67  : foundInf(false)
68  , foundNan(false)
69  , foundZeroDiag(false)
70  , foundZeroRowNorm(false) {}
71 
72  EquilibrationInfo(const std::size_t lclNumRows,
73  const std::size_t lclNumCols,
74  const bool assumeSymmetric_)
75  : rowNorms(Kokkos::View<mag_type*, device_type>("rowNorms", lclNumRows))
76  , rowDiagonalEntries(Kokkos::View<val_type*, device_type>("rowDiagonalEntries", lclNumRows))
77  , colNorms(Kokkos::View<mag_type*, device_type>("colNorms", lclNumCols))
78  , colDiagonalEntries(Kokkos::View<val_type*, device_type>("colDiagonalEntries",
79  assumeSymmetric_ ? std::size_t(0) : lclNumCols))
80  , rowScaledColNorms(Kokkos::View<mag_type*, device_type>("rowScaledColNorms",
81  assumeSymmetric_ ? std::size_t(0) : lclNumCols))
82  , assumeSymmetric(assumeSymmetric_)
83  , foundInf(false)
84  , foundNan(false)
85  , foundZeroDiag(false)
86  , foundZeroRowNorm(false) {}
87 
88  EquilibrationInfo(const Kokkos::View<mag_type*, device_type>& rowNorms_,
89  const Kokkos::View<val_type*, device_type>& rowDiagonalEntries_,
90  const Kokkos::View<mag_type*, device_type>& colNorms_,
91  const Kokkos::View<val_type*, device_type>& colDiagonalEntries_,
92  const Kokkos::View<mag_type*, device_type>& rowScaledColNorms_,
93  const bool assumeSymmetric_,
94  const bool foundInf_,
95  const bool foundNan_,
96  const bool foundZeroDiag_,
97  const bool foundZeroRowNorm_)
98  : rowNorms(rowNorms_)
99  , rowDiagonalEntries(rowDiagonalEntries_)
100  , colNorms(colNorms_)
101  , colDiagonalEntries(colDiagonalEntries_)
102  , rowScaledColNorms(rowScaledColNorms_)
103  , assumeSymmetric(assumeSymmetric_)
104  , foundInf(foundInf_)
105  , foundNan(foundNan_)
106  , foundZeroDiag(foundZeroDiag_)
107  , foundZeroRowNorm(foundZeroRowNorm_) {}
108 
110  template <class SrcDeviceType>
111  void
113  using execution_space = typename device_type::execution_space;
114  // DEEP_COPY REVIEW - DEVICE-TO-DEVICE
115  Kokkos::deep_copy(execution_space(), rowNorms, src.rowNorms);
116  // DEEP_COPY REVIEW - DEVICE-TO-DEVICE
118  // DEEP_COPY REVIEW - DEVICE-TO-DEVICE
119  Kokkos::deep_copy(execution_space(), colNorms, src.colNorms);
120  if (src.colDiagonalEntries.extent(0) == 0) {
122  Kokkos::View<val_type*, device_type>("colDiagonalEntries", 0);
123  } else {
124  // DEEP_COPY REVIEW - DEVICE-TO-DEVICE
126  }
127  if (src.rowScaledColNorms.extent(0) == 0) {
129  Kokkos::View<mag_type*, device_type>("rowScaledColNorms", 0);
130  } else {
131  // DEEP_COPY REVIEW - DEVICE-TO-DEVICE
132  Kokkos::deep_copy(execution_space(), rowScaledColNorms, src.rowScaledColNorms);
133  }
134 
136  foundInf = src.foundInf;
137  foundNan = src.foundNan;
140  }
141 
143  createMirrorView() {
144  auto rowNorms_h = Kokkos::create_mirror_view(rowNorms);
145  auto rowDiagonalEntries_h = Kokkos::create_mirror_view(rowDiagonalEntries);
146  auto colNorms_h = Kokkos::create_mirror_view(colNorms);
147  auto colDiagonalEntries_h = Kokkos::create_mirror_view(colDiagonalEntries);
148  auto rowScaledColNorms_h = Kokkos::create_mirror_view(rowScaledColNorms);
149 
150  return HostMirror{rowNorms_h, rowDiagonalEntries_h, colNorms_h,
151  colDiagonalEntries_h, rowScaledColNorms_h, assumeSymmetric,
153  }
154 
155  // We call a row a "diagonally dominant row" if the absolute value
156  // of the diagonal entry is >= the sum of the absolute values of the
157  // off-diagonal entries. The row norm is the sum of those two
158  // things, so this means diagAbsVal >= rowNorm - diagAbsVal. Ditto
159  // for a column.
160 
162  Kokkos::View<mag_type*, device_type> rowNorms;
163 
165  Kokkos::View<val_type*, device_type> rowDiagonalEntries;
166 
172  Kokkos::View<mag_type*, device_type> colNorms;
173 
177  Kokkos::View<val_type*, device_type> colDiagonalEntries;
178 
189  Kokkos::View<mag_type*, device_type> rowScaledColNorms;
190 
196 
198  bool foundInf;
199 
201  bool foundNan;
202 
205 
208 };
209 
210 } // namespace Details
211 } // namespace Tpetra
212 
213 #endif // TPETRA_DETAILS_EQUILIBRATIONINFO_HPP
Kokkos::View< mag_type *, device_type > colNorms
One-norms of the matrix&#39;s columns, distributed via the column Map.
Kokkos::View< mag_type *, device_type > rowNorms
One-norms of the matrix&#39;s rows, distributed via the row Map.
bool foundZeroDiag
Found a zero diagonal entry somewhere in the matrix.
bool assumeSymmetric
Whether to assume that the matrix is (globally) symmetric.
Kokkos::View< val_type *, device_type > colDiagonalEntries
Diagonal entries of the matrix, distributed via the column Map.
Struct storing results of Tpetra::computeRowAndColumnOneNorms.
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst.
Kokkos::View< val_type *, device_type > rowDiagonalEntries
Diagonal entries of the matrix, distributed via the row Map.
Kokkos::View< mag_type *, device_type > rowScaledColNorms
One-norms of the matrix&#39;s columns, after the matrix&#39;s rows have been scaled by rowNorms.
void assign(const EquilibrationInfo< ScalarType, SrcDeviceType > &src)
Deep-copy src into *this.
bool foundNan
Found a NaN somewhere in the matrix.
bool foundInf
Found an Inf somewhere in the matrix.
bool foundZeroRowNorm
At least one row of the matrix has a zero norm.