10#ifndef TPETRA_DETAILS_NORMIMPL_HPP
11#define TPETRA_DETAILS_NORMIMPL_HPP
22#include "TpetraCore_config.h"
23#include "Kokkos_Core.hpp"
24#include "Teuchos_ArrayView.hpp"
25#include "Teuchos_CommHelpers.hpp"
26#include "KokkosBlas.hpp"
27#include "KokkosKernels_ArithTraits.hpp"
30#ifndef DOXYGEN_SHOULD_SKIP_THIS
34template <
class OrdinalType>
59 const Kokkos::View<const ValueType**, ArrayLayout, DeviceType>&
X,
61 const Teuchos::ArrayView<const size_t>&
whichVecs,
62 const bool isConstantStride,
63 const bool isDistributed,
64 const Teuchos::Comm<int>* comm);
77template <
class RV,
class XMV>
78void lclNormImpl(
const RV& normsOut,
81 const Teuchos::ArrayView<const size_t>& whichVecs,
82 const bool constantStride,
85 using Kokkos::subview;
86 using mag_type =
typename RV::non_const_value_type;
88 Details::ProfilingRegion region(
"Tpetra::Details::lclNormImpl");
90 static_assert(
static_cast<int>(RV::rank) == 1,
91 "Tpetra::MultiVector::lclNormImpl: "
92 "The first argument normsOut must have rank 1.");
93 static_assert(Kokkos::is_view<XMV>::value,
94 "Tpetra::MultiVector::lclNormImpl: "
95 "The second argument X is not a Kokkos::View.");
96 static_assert(
static_cast<int>(XMV::rank) == 2,
97 "Tpetra::MultiVector::lclNormImpl: "
98 "The second argument X must have rank 2.");
100 const size_t lclNumRows =
static_cast<size_t>(X.extent(0));
101 TEUCHOS_TEST_FOR_EXCEPTION(lclNumRows != 0 && constantStride &&
102 static_cast<size_t>(X.extent(1)) != numVecs,
103 std::logic_error,
"Constant Stride X's dimensions are " << X.extent(0) <<
" x " << X.extent(1) <<
", which differ from the local dimensions " << lclNumRows <<
" x " << numVecs <<
". Please report this bug to "
104 "the Tpetra developers.");
105 TEUCHOS_TEST_FOR_EXCEPTION(lclNumRows != 0 && !constantStride &&
106 static_cast<size_t>(X.extent(1)) < numVecs,
107 std::logic_error,
"Strided X's dimensions are " << X.extent(0) <<
" x " << X.extent(1) <<
", which are incompatible with the local dimensions " << lclNumRows <<
" x " << numVecs <<
". Please report this bug to "
108 "the Tpetra developers.");
110 if (lclNumRows == 0) {
111 const mag_type zeroMag = KokkosKernels::ArithTraits<mag_type>::zero();
113 using execution_space =
typename RV::execution_space;
114 Kokkos::deep_copy(execution_space(), normsOut, zeroMag);
116 if (constantStride) {
117 if (whichNorm == NORM_INF) {
118 KokkosBlas::nrminf(normsOut, X);
119 }
else if (whichNorm == NORM_ONE) {
120 KokkosBlas::nrm1(normsOut, X);
121 }
else if (whichNorm == NORM_TWO) {
122 KokkosBlas::nrm2_squared(normsOut, X);
124 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Should never get here!");
132 for (
size_t k = 0; k < numVecs; ++k) {
133 const size_t X_col = constantStride ? k : whichVecs[k];
134 if (whichNorm == NORM_INF) {
135 KokkosBlas::nrminf(subview(normsOut, k),
136 subview(X, ALL(), X_col));
137 }
else if (whichNorm == NORM_ONE) {
138 KokkosBlas::nrm1(subview(normsOut, k),
139 subview(X, ALL(), X_col));
140 }
else if (whichNorm == NORM_TWO) {
141 KokkosBlas::nrm2_squared(subview(normsOut, k),
142 subview(X, ALL(), X_col));
144 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Should never get here!");
153template <
class ViewType>
154class SquareRootFunctor {
156 typedef typename ViewType::execution_space execution_space;
157 typedef typename ViewType::size_type size_type;
159 SquareRootFunctor(
const ViewType& theView)
160 : theView_(theView) {}
162 KOKKOS_INLINE_FUNCTION
void
163 operator()(
const size_type& i)
const {
164 typedef typename ViewType::non_const_value_type value_type;
165 typedef KokkosKernels::ArithTraits<value_type> KAT;
166 theView_(i) = KAT::sqrt(theView_(i));
174void gblNormImpl(
const RV& normsOut,
175 const Teuchos::Comm<int>*
const comm,
176 const bool distributed,
178 using Teuchos::REDUCE_MAX;
179 using Teuchos::REDUCE_SUM;
180 using Teuchos::reduceAll;
181 typedef typename RV::non_const_value_type mag_type;
183 Details::ProfilingRegion region(
"Tpetra::Details::gblNormImpl");
185 const size_t numVecs = normsOut.extent(0);
200 if (distributed && comm !=
nullptr) {
204 const int nv =
static_cast<int>(numVecs);
207 if (commIsInterComm) {
208 RV lclNorms(Kokkos::ViewAllocateWithoutInitializing(
"MV::normImpl lcl"), numVecs);
210 using execution_space =
typename RV::execution_space;
211 Kokkos::deep_copy(execution_space(), lclNorms, normsOut);
212 const mag_type*
const lclSum = lclNorms.data();
213 mag_type*
const gblSum = normsOut.data();
215 if (whichNorm == NORM_INF) {
216 reduceAll<int, mag_type>(*comm, REDUCE_MAX, nv, lclSum, gblSum);
218 reduceAll<int, mag_type>(*comm, REDUCE_SUM, nv, lclSum, gblSum);
221 mag_type*
const gblSum = normsOut.data();
222 if (whichNorm == NORM_INF) {
223 reduceAll<int, mag_type>(*comm, REDUCE_MAX, nv, gblSum, gblSum);
225 reduceAll<int, mag_type>(*comm, REDUCE_SUM, nv, gblSum, gblSum);
230 if (whichNorm == NORM_TWO) {
236 const bool inHostMemory =
237 std::is_same<
typename RV::memory_space,
238 typename RV::host_mirror_space::memory_space>::value;
240 for (
size_t j = 0; j < numVecs; ++j) {
241 normsOut(j) = KokkosKernels::ArithTraits<mag_type>::sqrt(normsOut(j));
248 SquareRootFunctor<RV> f(normsOut);
249 typedef typename RV::execution_space execution_space;
250 typedef Kokkos::RangePolicy<execution_space, size_t> range_type;
251 Kokkos::parallel_for(range_type(0, numVecs), f);
258template <
class ValueType,
263 const Kokkos::View<const ValueType**, ArrayLayout, DeviceType>&
X,
265 const Teuchos::ArrayView<const size_t>&
whichVecs,
266 const bool isConstantStride,
267 const bool isDistributed,
268 const Teuchos::Comm<int>* comm) {
269 using execution_space =
typename DeviceType::execution_space;
270 using RV = Kokkos::View<MagnitudeType*, Kokkos::HostSpace>;
274 const size_t numVecs = isConstantStride ?
static_cast<size_t>(
X.extent(1)) :
static_cast<size_t>(
whichVecs.size());
Declaration of Tpetra::Details::Profiling, a scope guard for Kokkos Profiling.
Struct that holds views of the contents of a CrsMatrix.
Implementation details of Tpetra.
void normImpl(MagnitudeType norms[], const Kokkos::View< const ValueType **, ArrayLayout, DeviceType > &X, const EWhichNorm whichNorm, const Teuchos::ArrayView< const size_t > &whichVecs, const bool isConstantStride, const bool isDistributed, const Teuchos::Comm< int > *comm)
Implementation of MultiVector norms.
EWhichNorm
Input argument for normImpl() (which see).
bool isInterComm(const Teuchos::Comm< int > &)
Return true if and only if the input communicator wraps an MPI intercommunicator.
Namespace Tpetra contains the class and methods constituting the Tpetra library.