Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_BlockedTpetraOperator.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_BlockedTpetraOperator_hpp__
11 #define __Teko_BlockedTpetraOperator_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_TpetraBlockedMappingStrategy.hpp"
25 #include "Teko_ConfigDefs.hpp"
26 
27 namespace Teko {
28 namespace TpetraHelpers {
29 
34  public:
43  BlockedTpetraOperator(const std::vector<std::vector<GO> > &vars,
44  const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > &content,
45  const std::string &label = "<ANYM>");
46 
55  virtual void SetContent(const std::vector<std::vector<GO> > &vars,
56  const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > &content);
57 
61  virtual void RebuildOps() { BuildBlockedOperator(); }
62 
63  virtual const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > GetContent() const {
64  return fullContent_;
65  }
66 
67  virtual const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > GetContent() {
68  return fullContent_;
69  }
70 
71  const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > GetBlock(int i, int j) const;
72 
76  void Reorder(const BlockReorderManager &brm);
77 
79  void RemoveReording();
80 
83  virtual void WriteBlocks(const std::string &prefix) const;
84 
85  // functions overloading Tpetra::Operator<ST,LO,GO,NT>
87 
88  // destructor
89  virtual ~BlockedTpetraOperator() {}
90 
91  // attribute set methods
92 
93  // don't use transpose...ever!
94  virtual int SetUseTranspose(bool /* useTranspose */) { return -1; }
95 
96  virtual int ApplyInverse(const Tpetra::MultiVector<ST, LO, GO, NT> & /* X */,
97  Tpetra::MultiVector<ST, LO, GO, NT> & /* Y */) const {
98  TEUCHOS_ASSERT(false);
99  return -1;
100  }
101 
102  virtual double NormInf() const {
103  TEUCHOS_ASSERT(false);
104  return 0.0;
105  }
106 
107  // attribute access functions
108  virtual bool UseTranspose() const { return false; }
109  virtual bool HasNormInf() const { return false; }
110  virtual const Teuchos::Comm<int> &Comm() const { return *fullContent_->getRangeMap()->getComm(); }
111 
113  bool testAgainstFullOperator(int count, ST tol) const;
114 
115  protected:
116  // gooey center of this shell
117  Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > fullContent_;
118  Teuchos::RCP<TpetraBlockedMappingStrategy> blockedMapping_;
119  Teuchos::RCP<Thyra::LinearOpBase<ST> > blockedOperator_;
120  Teuchos::RCP<const BlockReorderManager> reorderManager_;
121 
122  std::string label_;
123 
124  void BuildBlockedOperator();
125 };
126 
127 } // end namespace TpetraHelpers
128 } // end namespace Teko
129 
130 #endif
void Reorder(const BlockReorderManager &brm)
void RemoveReording()
Remove any reordering on this object.
Tear about a user specified Tpetra::Operator&lt;ST,LO,GO,NT&gt; (CrsMatrix) using a vector of vectors of GI...
virtual void WriteBlocks(const std::string &prefix) const
BlockedTpetraOperator(const std::vector< std::vector< GO > > &vars, const Teuchos::RCP< const Tpetra::Operator< ST, LO, GO, NT > > &content, const std::string &label="<ANYM>")
bool testAgainstFullOperator(int count, ST tol) const
Helps perform sanity checks.
virtual void SetContent(const std::vector< std::vector< GO > > &vars, const Teuchos::RCP< const Tpetra::Operator< ST, LO, GO, NT > > &content)
Implements the Epetra_Operator interface with a Thyra LinearOperator. This enables the use of absrtac...