Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_LinearProblem_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_LINEARPROBLEM_DECL_HPP
11#define TPETRA_LINEARPROBLEM_DECL_HPP
12
15
16#include "Teuchos_DataAccess.hpp"
17
20#include "Tpetra_RowMatrix_decl.hpp"
21#include "Tpetra_DistObject.hpp"
22#include "Tpetra_Details_ExecutionSpacesUser.hpp"
23
24namespace Tpetra {
25
42
43template <class Scalar,
44 class LocalOrdinal,
45 class GlobalOrdinal,
46 class Node>
47class LinearProblem : public DistObject<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
49 private:
52
53 public:
55
56
62
64
66
67
75
79 LinearProblem(const Teuchos::RCP<row_matrix_type>& A,
80 const Teuchos::RCP<multivector_type>& X,
81 const Teuchos::RCP<multivector_type>& B);
82
85 GlobalOrdinal, Node>& Problem);
86
88 virtual ~LinearProblem() = default;
89
91
93
94
101 void checkInput() const;
102
104
106
107
108 virtual bool
109 checkSizes(const SrcDistObject& source) override;
110
112
114
115
119 void setMatrix(Teuchos::RCP<row_matrix_type> A) { A_ = A; }
120
124 void setLHS(Teuchos::RCP<multivector_type> X) { X_ = X; }
125
129 void setRHS(Teuchos::RCP<multivector_type> B) { B_ = B; }
130
132
134
135
147 void leftScale(const Teuchos::RCP<const vector_type>& D,
148 Teuchos::ETransp mode = Teuchos::NO_TRANS);
149
163 void rightScale(const Teuchos::RCP<const vector_type>& D,
164 Teuchos::ETransp mode = Teuchos::NO_TRANS);
165
167
169
170
172 Teuchos::RCP<row_matrix_type> getMatrix() const { return (A_); }
174 Teuchos::RCP<multivector_type> getLHS() const { return (X_); }
176 Teuchos::RCP<multivector_type> getRHS() const { return (B_); }
177
179
180 private:
181 Teuchos::RCP<row_matrix_type> A_;
182 Teuchos::RCP<multivector_type> X_;
183 Teuchos::RCP<multivector_type> B_;
184
186 GlobalOrdinal, Node>& Problem) = default;
187};
188
189} // namespace Tpetra
190
191#endif // TPETRA_LINEARPROBLEM_DECL_HPP
Declaration of the Tpetra::MultiVector class.
Declaration of the Tpetra::Vector class.
Struct that holds views of the contents of a CrsMatrix.
A class can inherit from this if it wants to use Tpetra managed spaces.
Base class for distributed Tpetra objects that support data redistribution.
Class that encapulates linear problem (Ax = b).
virtual bool checkSizes(const SrcDistObject &source) override
Compare the source and target (this) objects for compatibility.
Teuchos::RCP< row_matrix_type > getMatrix() const
Get an RCP to the matrix A.
void rightScale(const Teuchos::RCP< const vector_type > &D, Teuchos::ETransp mode=Teuchos::NO_TRANS)
Perform right scaling of a linear problem.
Teuchos::RCP< multivector_type > getLHS() const
Get an RCP to the left-hand-side X.
LinearProblem()
Default Constructor.
Teuchos::RCP< multivector_type > getRHS() const
Get an RCP to the right-hand-side B.
void setMatrix(Teuchos::RCP< row_matrix_type > A)
Set Matrix A of linear problem AX = B using a RowMatrix.
void leftScale(const Teuchos::RCP< const vector_type > &D, Teuchos::ETransp mode=Teuchos::NO_TRANS)
Perform left scaling of a linear problem.
void setRHS(Teuchos::RCP< multivector_type > B)
Set right-hand-side B of linear problem AX = B.
void checkInput() const
Check input parameters for existence and size consistency.
void setLHS(Teuchos::RCP< multivector_type > X)
Set left-hand-side X of linear problem AX = B.
virtual ~LinearProblem()=default
LinearProblem Destructor.
One or more distributed dense vectors.
A read-only, row-oriented interface to a sparse matrix.
Abstract base class for objects that can be the source of an Import or Export operation.
A distributed dense vector.
Namespace Tpetra contains the class and methods constituting the Tpetra library.