Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Details_scaleBlockDiagonal.hpp
Go to the documentation of this file.
1// @HEADER
2// ***********************************************************************
3//
4// Tpetra: Templated Linear Algebra Services Package
5// Copyright (2008) Sandia Corporation
6//
7// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8// the U.S. Government retains certain rights in this software.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38//
39// ************************************************************************
40// @HEADER
41
42#ifndef TPETRA_DETAILS_SCALEBLOCKDIAGONAL_HPP
43#define TPETRA_DETAILS_SCALEBLOCKDIAGONAL_HPP
44
45#include "TpetraCore_config.h"
46#include "Tpetra_CrsMatrix.hpp"
47#include "Teuchos_ScalarTraits.hpp"
49#include "KokkosBatched_Util.hpp"
50#include "KokkosBatched_LU_Decl.hpp"
51#include "KokkosBatched_LU_Serial_Impl.hpp"
52#include "KokkosBatched_Trsm_Decl.hpp"
53#include "KokkosBatched_Trsm_Serial_Impl.hpp"
54
55
62
63namespace Tpetra {
64namespace Details {
65
66
67template<class MultiVectorType>
68void inverseScaleBlockDiagonal(MultiVectorType & blockDiagonal, bool doTranspose, MultiVectorType & multiVectorToBeScaled) {
69 using LO = typename MultiVectorType::local_ordinal_type;
70 using range_type = Kokkos::RangePolicy<typename MultiVectorType::node_type::execution_space, LO>;
71 using namespace KokkosBatched;
72 typename MultiVectorType::impl_scalar_type SC_one = Teuchos::ScalarTraits<typename MultiVectorType::impl_scalar_type>::one();
73
74 // Sanity checking: Map Compatibility (A's rowmap matches diagonal's map)
76 TEUCHOS_TEST_FOR_EXCEPTION(!blockDiagonal.getMap()->isSameAs(*multiVectorToBeScaled.getMap()),
77 std::runtime_error, "Tpetra::Details::scaledBlockDiagonal was given incompatible maps");
78 }
79
80 LO numrows = blockDiagonal.getLocalLength();
81 LO blocksize = blockDiagonal.getNumVectors();
82 LO numblocks = numrows / blocksize;
83
84 // Get Kokkos versions of objects
85
86 auto blockDiag = blockDiagonal.getLocalViewDevice(Access::OverwriteAll);
87 auto toScale = multiVectorToBeScaled.getLocalViewDevice(Access::ReadWrite);
88
89 typedef Algo::Level3::Unblocked algo_type;
90 Kokkos::parallel_for("scaleBlockDiagonal",range_type(0,numblocks),KOKKOS_LAMBDA(const LO i){
91 Kokkos::pair<LO,LO> row_range(i*blocksize,(i+1)*blocksize);
92 auto A = Kokkos::subview(blockDiag,row_range, Kokkos::ALL());
93 auto B = Kokkos::subview(toScale, row_range, Kokkos::ALL());
94
95 // Factor
96 SerialLU<algo_type>::invoke(A);
97
98 if(doTranspose) {
99 // Solve U^T
100 SerialTrsm<Side::Left,Uplo::Upper,Trans::Transpose,Diag::NonUnit,algo_type>::invoke(SC_one,A,B);
101 // Solver L^T
102 SerialTrsm<Side::Left,Uplo::Lower,Trans::Transpose,Diag::Unit,algo_type>::invoke(SC_one,A,B);
103 }
104 else {
105 // Solve L
106 SerialTrsm<Side::Left,Uplo::Lower,Trans::NoTranspose,Diag::Unit,algo_type>::invoke(SC_one,A,B);
107 // Solve U
108 SerialTrsm<Side::Left,Uplo::Upper,Trans::NoTranspose,Diag::NonUnit,algo_type>::invoke(SC_one,A,B);
109 }
110 });
111}
112
113} // namespace Details
114} // namespace Tpetra
115
116#endif // TPETRA_DETAILS_SCALEBLOCKDIAGONAL_HPP
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.
static bool debug()
Whether Tpetra is in debug mode.
Implementation details of Tpetra.
Namespace Tpetra contains the class and methods constituting the Tpetra library.