MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_BlockedRAPFactory_decl.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// MueLu: A package for multigrid based preconditioning
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact
39// Jonathan Hu (jhu@sandia.gov)
40// Andrey Prokopenko (aprokop@sandia.gov)
41// Ray Tuminaro (rstumin@sandia.gov)
42//
43// ***********************************************************************
44//
45// @HEADER
46#ifndef MUELU_BLOCKEDRAPFACTORY_DECL_HPP
47#define MUELU_BLOCKEDRAPFACTORY_DECL_HPP
48
49#include <Xpetra_Matrix_fwd.hpp>
50#include <Xpetra_CrsMatrix_fwd.hpp>
51#include <Xpetra_CrsMatrixWrap_fwd.hpp>
52#include <Xpetra_BlockedCrsMatrix_fwd.hpp>
53#include <Xpetra_MatrixFactory_fwd.hpp>
54
55#include "MueLu_ConfigDefs.hpp"
57
58#include "MueLu_Level_fwd.hpp"
63
64namespace MueLu {
69 template <class Scalar = DefaultScalar,
72 class Node = DefaultNode>
74#undef MUELU_BLOCKEDRAPFACTORY_SHORT
76
77 public:
79
80
82
83 virtual ~BlockedRAPFactory() = default;
85
87
88
89 RCP<const ParameterList> GetValidParameterList() const override;
90
91 void DeclareInput(Level &fineLevel, Level &coarseLevel) const override;
92
94
96
97 void Build(Level &fineLevel, Level &coarseLevel) const override;
99
101
102
104 void SetRepairZeroDiagonal(bool const &repair) {
105 repairZeroDiagonals_ = repair;
106 if(repair) checkAc_ = true; // make sure that plausibility check is performed. Otherwise SetRepairZeroDiagonal(true) has no effect.
107 }
108
110 void SetPlausibilityCheck(bool const &check) {
111 checkAc_ = check;
112 }
113
115
120 void AddTransferFactory(const RCP<const FactoryBase>& factory);
121
122 // TODO add a function to remove a specific transfer factory?
123
125 size_t NumTransferFactories() const { return transferFacts_.size(); }
126
128
129 private:
130
133 static void CheckMainDiagonal(RCP<BlockedCrsMatrix> & bAc, bool repairZeroDiagonals = false);
134
138
143
145
147 std::vector<RCP<const FactoryBase> > transferFacts_;
148
150
151 }; //class BlockedRAPFactory
152
153} //namespace MueLu
154
155#define MUELU_BLOCKEDRAPFACTORY_SHORT
156#endif // MUELU_BLOCKEDRAPFACTORY_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Factory for building coarse matrices.
void SetPlausibilityCheck(bool const &check)
Indicate that a simple plausibility check shall be done for Ac after building RAP.
void SetRepairZeroDiagonal(bool const &repair)
Indicate that zero entries on the diagonal of Ac shall be repaired (i.e. if A(i,i) == 0....
size_t NumTransferFactories() const
Returns number of transfer factories.
static void CheckMainDiagonal(RCP< BlockedCrsMatrix > &bAc, bool repairZeroDiagonals=false)
virtual ~BlockedRAPFactory()=default
void Build(Level &fineLevel, Level &coarseLevel) const override
Build an object with this factory.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const override
Input.
std::vector< RCP< const FactoryBase > > transferFacts_
list of user-defined transfer Factories
RCP< const ParameterList > GetValidParameterList() const override
Return a const parameter list of valid parameters that setParameterList() will accept.
void AddTransferFactory(const RCP< const FactoryBase > &factory)
Add transfer factory in the end of list of transfer factories in RepartitionAcFactory.
Class that holds all level-specific information.
Base class for factories that use two levels (fineLevel and coarseLevel).
Namespace for MueLu classes and methods.
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar