Teko Version of the Day
Loading...
Searching...
No Matches
Teko_DiagonalPreconditionerFactory.cpp
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#include "Teko_DiagonalPreconditionerFactory.hpp"
11#include "Teko_DiagonalPreconditionerOp.hpp"
12#ifdef TEKO_HAVE_EPETRA
13#include "Thyra_get_Epetra_Operator.hpp"
14#include "Epetra_CrsMatrix.h"
15#include "EpetraExt_PointToBlockDiagPermute.h"
16#endif
17
18#include "Teko_TpetraHelpers.hpp"
19#ifdef TEKO_HAVE_EPETRA
20#include "Thyra_EpetraLinearOp.hpp"
21#endif
22#include "Thyra_TpetraLinearOp.hpp"
23
24using Teuchos::rcp;
25using Teuchos::RCP;
26
27namespace Teko {
28
29DiagonalPrecondState::DiagonalPrecondState() {}
30
31/*****************************************************/
32
33DiagonalPreconditionerFactory::DiagonalPreconditionerFactory() {}
34
36 return Teuchos::rcp(new DiagonalPrecondState());
37}
38
40 LinearOp& lo, PreconditionerState& state) const {
41 if (diagonalType_ == BlkDiag) {
42 TEUCHOS_TEST_FOR_EXCEPTION(TpetraHelpers::isTpetraLinearOp(lo), std::runtime_error,
43 "BlkDiag not implemented for Tpetra operators");
44#ifdef TEKO_HAVE_EPETRA
45 // Sanity check the state
46 DiagonalPrecondState& MyState = Teuchos::dyn_cast<DiagonalPrecondState>(state);
47
48 // Get the underlying Epetra_CrsMatrix, if we have one
49 Teuchos::RCP<const Epetra_Operator> eo = Thyra::get_Epetra_Operator(*lo);
50 TEUCHOS_ASSERT(eo != Teuchos::null);
51 Teuchos::RCP<const Epetra_CrsMatrix> MAT =
52 Teuchos::rcp_dynamic_cast<const Epetra_CrsMatrix>(eo);
53 TEUCHOS_ASSERT(MAT != Teuchos::null);
54
55 // Create a new EpetraExt_PointToBlockDiagPermute for the state object, if we don't have one
56 Teuchos::RCP<EpetraExt_PointToBlockDiagPermute> BDP;
57 if (MyState.BDP_ == Teuchos::null) {
58 BDP = Teuchos::rcp(new EpetraExt_PointToBlockDiagPermute(*MAT));
59 BDP->SetParameters(List_);
60 BDP->Compute();
61 MyState.BDP_ = BDP;
62 }
63
64 RCP<Epetra_FECrsMatrix> Hcrs = rcp(MyState.BDP_->CreateFECrsMatrix());
65 return Thyra::epetraLinearOp(Hcrs);
66
67 // Build the LinearOp object (NTS: swapping the range and domain)
68 // LinearOp MyOp = Teuchos::rcp(new
69 // DiagonalPreconditionerOp(MyState.BDP_,lo->domain(),lo->range()));
70#endif
71 }
72
73 return getInvDiagonalOp(lo, diagonalType_);
74}
75
76void DiagonalPreconditionerFactory::initializeFromParameterList(const Teuchos::ParameterList& pl) {
77 List_ = pl;
78
79 diagonalType_ = BlkDiag;
80 if (pl.isParameter("Diagonal Type")) {
81 diagonalType_ = getDiagonalType(pl.get<std::string>("Diagonal Type"));
82 TEUCHOS_TEST_FOR_EXCEPT(diagonalType_ == NotDiag);
83 }
84
85 if (diagonalType_ == BlkDiag) {
86 // Reset default to invert mode if the user hasn't specified something else
87 Teuchos::ParameterList& SubList = List_.sublist("blockdiagmatrix: list");
88 SubList.set("apply mode", SubList.get("apply mode", "invert"));
89 }
90}
91
92} // end namespace Teko
@ BlkDiag
Specifies that a block diagonal approximation is to be used.
@ NotDiag
For user convenience, if Teko recieves this value, exceptions will be thrown.
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl)
Initialize from a parameter list.
Teuchos::RCP< PreconditionerState > buildPreconditionerState() const
Builds a preconditioner state object.
LinearOp buildPreconditionerOperator(LinearOp &lo, PreconditionerState &state) const
An implementation of a state object preconditioners.