Tpetra parallel linear algebra
Version of the Day
Toggle main menu visibility
Loading...
Searching...
No Matches
core
src
Tpetra_Details_getGraphOffRankOffsets_def.hpp
Go to the documentation of this file.
1
/*
2
// @HEADER
3
// ***********************************************************************
4
//
5
// Tpetra: Templated Linear Algebra Services Package
6
// Copyright (2008) Sandia Corporation
7
//
8
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9
// the U.S. Government retains certain rights in this software.
10
//
11
// Redistribution and use in source and binary forms, with or without
12
// modification, are permitted provided that the following conditions are
13
// met:
14
//
15
// 1. Redistributions of source code must retain the above copyright
16
// notice, this list of conditions and the following disclaimer.
17
//
18
// 2. Redistributions in binary form must reproduce the above copyright
19
// notice, this list of conditions and the following disclaimer in the
20
// documentation and/or other materials provided with the distribution.
21
//
22
// 3. Neither the name of the Corporation nor the names of the
23
// contributors may be used to endorse or promote products derived from
24
// this software without specific prior written permission.
25
//
26
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37
//
38
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39
//
40
// ************************************************************************
41
// @HEADER
42
*/
43
44
#ifndef TPETRA_DETAILS_GETGRAPHOFFRANKOFFSETS_DEF_HPP
45
#define TPETRA_DETAILS_GETGRAPHOFFRANKOFFSETS_DEF_HPP
46
51
52
#include "
Tpetra_Details_OrdinalTraits.hpp
"
53
#include "Tpetra_Map.hpp"
54
#include "KokkosSparse_findRelOffset.hpp"
55
56
namespace
Tpetra
{
57
namespace
Details
{
58
namespace
Impl {
59
73
template
<
class
LO,
74
class
GO,
75
class
Node,
76
class
OffsetType>
77
GetGraphOffRankOffsets<LO, GO, Node, OffsetType>::
78
GetGraphOffRankOffsets
(
const
offsets_type& OffRankOffsets,
79
const
local_map_type& lclColMap,
80
const
local_map_type& lclDomMap,
81
const
row_offsets_type& ptr,
82
const
lcl_col_inds_type& ind) :
83
OffRankOffsets_ (OffRankOffsets),
84
lclColMap_ (lclColMap),
85
lclDomMap_ (lclDomMap),
86
ptr_ (ptr),
87
ind_ (ind)
88
{
89
typedef
typename
device_type::execution_space
execution_space
;
90
typedef
Kokkos::RangePolicy<execution_space, LO> policy_type;
91
92
lclNumRows_ = ptr.extent(0)-1;
93
policy_type range (0, ptr.extent(0));
94
Kokkos::parallel_for (range, *
this
);
95
}
96
97
template
<
class
LO,
98
class
GO,
99
class
Node,
100
class
OffsetType>
101
KOKKOS_FUNCTION
void
102
GetGraphOffRankOffsets<LO, GO, Node, OffsetType>::
103
operator()
(
const
LO& lclRowInd)
const
104
{
105
const
LO INVALID =
106
Tpetra::Details::OrdinalTraits<LO>::invalid ();
107
108
if
(lclRowInd == lclNumRows_)
109
OffRankOffsets_[lclRowInd] = ptr_[lclRowInd];
110
else
{
111
// TODO: use parallel reduce
112
size_t
offset = ptr_[lclRowInd+1];
113
for
(
size_t
j = ptr_[lclRowInd]; j < ptr_[lclRowInd+1]; j++) {
114
const
LO lclColInd = ind_[j];
115
const
GO gblColInd = lclColMap_.getGlobalElement (lclColInd);
116
const
LO lclDomInd = lclDomMap_.getLocalElement (gblColInd);
117
if
((lclDomInd == INVALID) && (j < offset))
118
offset = j;
119
}
120
OffRankOffsets_[lclRowInd] = offset;
121
}
122
}
123
124
}
// namespace Impl
125
}
// namespace Details
126
}
// namespace Tpetra
127
128
// Explicit template instantiation macro for
129
// Tpetra::Details::Impl::GetGraphOffRankOffsets. NOT FOR USERS!!! Must
130
// be used inside the Tpetra namespace.
131
#define TPETRA_DETAILS_IMPL_GETGRAPHOFFRANKOFFSETS_INSTANT( LO, GO, NODE ) \
132
template class Details::Impl::GetGraphOffRankOffsets< LO, GO, NODE::device_type >;
133
134
#endif
// TPETRA_DETAILS_GETGRAPHOFFRANKOFFSETS_DEF_HPP
Tpetra_Details_OrdinalTraits.hpp
Import KokkosSparse::OrdinalTraits, a traits class for "invalid" (flag) values of integer types,...
Tpetra::execution_space
Tpetra::Details::Impl::GetGraphOffRankOffsets::GetGraphOffRankOffsets
GetGraphOffRankOffsets(const offsets_type &OffRankOffsets, const local_map_type &lclColMap, const local_map_type &lclDomMap, const row_offsets_type &ptr, const lcl_col_inds_type &ind)
Constructor; also runs the functor.
Definition
Tpetra_Details_getGraphOffRankOffsets_def.hpp:78
Tpetra::Details::Impl::GetGraphOffRankOffsets::operator()
KOKKOS_FUNCTION void operator()(const LO &lclRowInd) const
Kokkos::parallel_for loop body.
Definition
Tpetra_Details_getGraphOffRankOffsets_def.hpp:103
Tpetra::Details
Nonmember function that computes a residual Computes R = B - A * X.
Definition
Tpetra_KokkosRefactor_Details_MultiVectorLocalDeepCopy.hpp:51
Tpetra
Namespace Tpetra contains the class and methods constituting the Tpetra library.
Generated by
1.17.0