Teko Version of the Day
Loading...
Searching...
No Matches
Teko_StridedTpetraOperator.hpp
1// @HEADER
2// *****************************************************************************
3// Teko: A package for block and physics based preconditioning
4//
5// Copyright 2010 NTESS and the Teko contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef __Teko_StridedTpetraOperator_hpp__
11#define __Teko_StridedTpetraOperator_hpp__
12
13// Tpetra includes
14#include "Tpetra_Operator.hpp"
15
16// Teuchos includes
17#include "Teuchos_RCP.hpp"
18
19#include "Thyra_LinearOpBase.hpp"
20
21// Teko includes
22#include "Teko_BlockedReordering.hpp"
23#include "Teko_TpetraOperatorWrapper.hpp"
24#include "Teko_TpetraStridedMappingStrategy.hpp"
25#include "Teko_ConfigDefs.hpp"
26
27namespace Teko {
28namespace TpetraHelpers {
29
30class StridedTpetraOperator : public TpetraOperatorWrapper {
31 public:
32 enum eNormType { Inf, One, Frobenius };
33
34 StridedTpetraOperator(int numVars,
35 const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > &content,
36 const std::string &label = "<ANYM>");
37 StridedTpetraOperator(const std::vector<int> &vars,
38 const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > &content,
39 const std::string &label = "<ANYM>");
40
41 virtual void SetContent(const std::vector<int> &vars,
42 const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > &content);
43
44 virtual void RebuildOps() { BuildBlockedOperator(); }
45
46 virtual const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > GetContent() const {
47 return fullContent_;
48 }
49
50 // virtual const Teuchos::RCP<Tpetra::Operator<ST,LO,GO,NT> > GetContent()
51 // { return fullContent_; }
52
53 const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > GetBlock(int i, int j) const;
54
58 void Reorder(const BlockReorderManager &brm);
59
61 void RemoveReording();
62
65 virtual void WriteBlocks(const std::string &prefix) const;
66
74 virtual std::string PrintNorm(const eNormType &nrmType = Frobenius, const char newline = '\n');
75
76 // functions overloading Tpetra::Operator<ST,LO,GO,NT>
78
79 // destructor
80 virtual ~StridedTpetraOperator() {}
81
82 // attribute set methods
83
84 // don't use transpose...ever!
85 virtual int SetUseTranspose(bool /* useTranspose */) { return -1; }
86
87 virtual int applyInverse(const Tpetra::MultiVector<ST, LO, GO, NT> & /* X */,
88 Tpetra::MultiVector<ST, LO, GO, NT> & /* Y */) const {
89 TEUCHOS_ASSERT(false);
90 return -1;
91 }
92
93 virtual ST NormInf() const {
94 TEUCHOS_ASSERT(false);
95 return 0.0;
96 }
97
98 // attribute access functions
99 virtual bool UseTranspose() const { return false; }
100 virtual bool HasNormInf() const { return false; }
101 virtual const Teuchos::Comm<int> &Comm() const { return *fullContent_->getRangeMap()->getComm(); }
102
104 bool testAgainstFullOperator(int count, ST tol) const;
105
106 protected:
107 // gooey center of this shell
108 Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > fullContent_;
109 Teuchos::RCP<TpetraStridedMappingStrategy> stridedMapping_;
110 Teuchos::RCP<Thyra::LinearOpBase<ST> > stridedOperator_;
111 Teuchos::RCP<const BlockReorderManager> reorderManager_;
112
113 std::string label_;
114
115 void BuildBlockedOperator();
116};
117
118} // end namespace TpetraHelpers
119} // end namespace Teko
120
121#endif
Implements the Epetra_Operator interface with a Thyra LinearOperator. This enables the use of absrtac...