MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_PatternFactory_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // MueLu: A package for multigrid based preconditioning
4 //
5 // Copyright 2012 NTESS and the MueLu contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef MUELU_PATTERNFACTORY_DEF_HPP
11 #define MUELU_PATTERNFACTORY_DEF_HPP
12 
13 #include <Xpetra_Matrix.hpp>
14 #include <Xpetra_MatrixMatrix.hpp>
15 
17 
18 #include "MueLu_MasterList.hpp"
19 #include "MueLu_Monitor.hpp"
20 //#include "MueLu_Utilities.hpp"
21 
22 namespace MueLu {
23 
24 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
26  RCP<ParameterList> validParamList = rcp(new ParameterList());
27 
28 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
29  SET_VALID_ENTRY("emin: pattern order");
30 #undef SET_VALID_ENTRY
31 
32  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory for the matrix");
33  validParamList->set<RCP<const FactoryBase> >("P", Teuchos::null, "Generating factory for the matrix providing nonzero graph");
34 
35  return validParamList;
36 }
37 
38 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40  Input(coarseLevel, "P");
41 
42  const ParameterList& pL = GetParameterList();
43  if (pL.get<int>("emin: pattern order") > 0)
44  Input(fineLevel, "A");
45 }
46 
47 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
49  FactoryMonitor m(*this, "Ppattern", coarseLevel);
50 
51  RCP<Matrix> P = Get<RCP<Matrix> >(coarseLevel, "P");
52 
53  const ParameterList& pL = GetParameterList();
54  int k = pL.get<int>("emin: pattern order");
55 
56  if (k > 0) {
57  RCP<Matrix> A = Get<RCP<Matrix> >(fineLevel, "A");
58  RCP<Matrix> AP;
59 
60  bool doFillComplete = true;
61  bool optimizeStorage = true;
62 
63  for (int i = 0; i < k; i++) {
64  AP = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*A, false, *P, false, GetOStream(Statistics2), doFillComplete, optimizeStorage);
65  P.swap(AP);
66  }
67  }
68 
69  Set(coarseLevel, "Ppattern", P->getCrsGraph());
70 }
71 
72 } // namespace MueLu
73 
74 #endif // MUELU_PATTERNFACTORY_DEF_HPP
#define SET_VALID_ENTRY(name)
T & get(const std::string &name, T def_value)
Timer to be used in factories. Similar to Monitor but with additional timers.
void swap(RCP< T > &r_ptr)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
Print even more statistics.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
static void Multiply(const Matrix &A, bool transposeA, const Matrix &B, bool transposeB, Matrix &C, bool call_FillComplete_on_result=true, bool doOptimizeStorage=true, const std::string &label=std::string(), const RCP< ParameterList > &params=null)
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.