Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_computeRowAndColumnOneNorms_decl.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_COMPUTEROWANDCOLUMNONENORMS_DECL_HPP
11#define TPETRA_COMPUTEROWANDCOLUMNONENORMS_DECL_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
24
25namespace Tpetra {
26
39template <class SC, class LO, class GO, class NT>
40#if KOKKOS_VERSION >= 40799
41Details::EquilibrationInfo<typename KokkosKernels::ArithTraits<SC>::val_type,
42#else
43Details::EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
44#endif
45 typename NT::device_type>
47
74template <class SC, class LO, class GO, class NT>
75#if KOKKOS_VERSION >= 40799
76Details::EquilibrationInfo<typename KokkosKernels::ArithTraits<SC>::val_type,
77#else
78Details::EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
79#endif
80 typename NT::device_type>
82 const bool assumeSymmetric);
83
84} // namespace Tpetra
85
86#endif // TPETRA_COMPUTEROWANDCOLUMNONENORMS_DECL_HPP
Declaration of Tpetra::Details::EquilibrationInfo.
Forward declaration of Tpetra::RowMatrix.
Struct that holds views of the contents of a CrsMatrix.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
Details::EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeRowOneNorms(const Tpetra::RowMatrix< SC, LO, GO, NT > &A)
Compute global row one-norms ("row sums") of the input sparse matrix A, in a way suitable for one-sid...
Details::EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeRowAndColumnOneNorms(const Tpetra::RowMatrix< SC, LO, GO, NT > &A, const bool assumeSymmetric)
Compute global row and column one-norms ("row sums" and "column sums") of the input sparse matrix A,...