Teko Version of the Day
Loading...
Searching...
No Matches
Teko_BlockInvDiagonalStrategy.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_BlockInvDiagonalStrategy_hpp__
11#define __Teko_BlockInvDiagonalStrategy_hpp__
12
13#include <vector>
14
15// Teuchos includes
16#include "Teuchos_RCP.hpp"
17
18// Thyra includes
19#include "Thyra_LinearOpBase.hpp"
20
21// Teko includes
22#include "Teko_Utilities.hpp"
23#include "Teko_InverseFactory.hpp"
24#include "Teko_BlockPreconditionerFactory.hpp"
25
26namespace Teko {
27
43 public:
44 virtual ~BlockInvDiagonalStrategy() {}
45
47 virtual void getInvD(const BlockedLinearOp &A, BlockPreconditionerState &state,
48 std::vector<LinearOp> &invDiag) const = 0;
49};
50
57class StaticInvDiagStrategy : public BlockInvDiagonalStrategy {
58 public:
59 StaticInvDiagStrategy(const LinearOp &invD0, const LinearOp &invD1) {
60 invDiag_.push_back(invD0);
61 invDiag_.push_back(invD1);
62 }
63
64 StaticInvDiagStrategy(const std::vector<LinearOp> &invD) : invDiag_(invD) {}
65
66 virtual ~StaticInvDiagStrategy() {}
67
71 virtual void getInvD(const BlockedLinearOp & /* A */, BlockPreconditionerState & /* state */,
72 std::vector<LinearOp> &invDiag) const {
73 invDiag.clear();
74 invDiag = invDiag_;
75 }
76
77 protected:
78 // stored inverse operators
79 std::vector<Teuchos::RCP<const Thyra::LinearOpBase<double> > > invDiag_;
80};
81
87 public:
93 InvFactoryDiagStrategy(const Teuchos::RCP<InverseFactory> &factory);
94
102 InvFactoryDiagStrategy(const std::vector<Teuchos::RCP<InverseFactory> > &factories,
103 const Teuchos::RCP<InverseFactory> &defaultFact = Teuchos::null);
104
106 const std::vector<Teuchos::RCP<InverseFactory> > &inverseFactories,
107 const std::vector<Teuchos::RCP<InverseFactory> > &preconditionerFactories,
108 const Teuchos::RCP<InverseFactory> &defaultInverseFact = Teuchos::null,
109 const Teuchos::RCP<InverseFactory> &defaultPreconditionerFact = Teuchos::null);
110
111 virtual ~InvFactoryDiagStrategy() {}
112
120 virtual void getInvD(const BlockedLinearOp &A, BlockPreconditionerState &state,
121 std::vector<LinearOp> &invDiag) const;
122
124 const std::vector<Teuchos::RCP<InverseFactory> > &getFactories() const { return invDiagFact_; }
125
126 protected:
127 // stored inverse operators
128 std::vector<Teuchos::RCP<InverseFactory> > invDiagFact_;
129 std::vector<Teuchos::RCP<InverseFactory> > precDiagFact_;
130 Teuchos::RCP<InverseFactory> defaultInvFact_;
131 Teuchos::RCP<InverseFactory> defaultPrecFact_;
132
134 LinearOp buildInverse(const InverseFactory &invFact, Teuchos::RCP<InverseFactory> &precFact,
135 const LinearOp &matrix, BlockPreconditionerState &state,
136 const std::string &opPrefix, int i) const;
137
138 private:
141};
142
143} // end namespace Teko
144
145#endif
virtual void getInvD(const BlockedLinearOp &A, BlockPreconditionerState &state, std::vector< LinearOp > &invDiag) const =0
returns an (approximate) inverse of the diagonal blocks of A
An implementation of a state object for block preconditioners.
InvFactoryDiagStrategy(const Teuchos::RCP< InverseFactory > &factory)
LinearOp buildInverse(const InverseFactory &invFact, Teuchos::RCP< InverseFactory > &precFact, const LinearOp &matrix, BlockPreconditionerState &state, const std::string &opPrefix, int i) const
Conveinence function for building inverse operators.
const std::vector< Teuchos::RCP< InverseFactory > > & getFactories() const
Get factories for testing purposes.
virtual void getInvD(const BlockedLinearOp &A, BlockPreconditionerState &state, std::vector< LinearOp > &invDiag) const
Abstract class for building an inverse operator.
virtual void getInvD(const BlockedLinearOp &, BlockPreconditionerState &, std::vector< LinearOp > &invDiag) const