MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_TentativePFactory_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_TENTATIVEPFACTORY_DECL_HPP
47#define MUELU_TENTATIVEPFACTORY_DECL_HPP
48
49#include <Teuchos_ScalarTraits.hpp>
50#include <Teuchos_SerialDenseMatrix.hpp>
51#include <Teuchos_SerialQRDenseSolver.hpp>
52
53#include <Xpetra_CrsGraphFactory_fwd.hpp>
54#include <Xpetra_CrsMatrix_fwd.hpp>
55#include <Xpetra_Matrix_fwd.hpp>
56#include <Xpetra_MultiVector_fwd.hpp>
57#include <Xpetra_MapFactory_fwd.hpp>
58#include <Xpetra_Map_fwd.hpp>
59#include <Xpetra_Matrix_fwd.hpp>
60#include <Xpetra_MultiVector_fwd.hpp>
61#include <Xpetra_MultiVectorFactory_fwd.hpp>
62#include <Xpetra_Import_fwd.hpp>
63#include <Xpetra_ImportFactory_fwd.hpp>
64#include <Xpetra_CrsMatrixWrap_fwd.hpp>
65
66#include "MueLu_ConfigDefs.hpp"
68
73#include "MueLu_Level_fwd.hpp"
75#include "MueLu_PFactory.hpp"
77
78namespace MueLu {
79
119template <class Scalar = DefaultScalar,
122 class Node = DefaultNode>
124#undef MUELU_TENTATIVEPFACTORY_SHORT
126
127 public:
129
130
133
135 virtual ~TentativePFactory() { }
137
138 RCP<const ParameterList> GetValidParameterList() const;
139
141
142
143 void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
144
146
148
149
150 void Build (Level& fineLevel, Level& coarseLevel) const;
151 void BuildP(Level& fineLevel, Level& coarseLevel) const;
152
154
155 private:
156
157 void BuildPuncoupled(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
158 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace, const int levelID) const;
159 void BuildPcoupled (RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
160 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace) const;
161 void BuildPuncoupledBlockCrs(RCP<Matrix> A, RCP<Aggregates> aggregates, RCP<AmalgamationInfo> amalgInfo, RCP<MultiVector> fineNullspace,
162 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative, RCP<MultiVector>& coarseNullspace, const int levelID) const;
163
164 mutable bool bTransferCoordinates_ = false;
165
166 }; //class TentativePFactory
167
168} //namespace MueLu
169
170//TODO: noQR_
171
172#define MUELU_TENTATIVEPFACTORY_SHORT
173#endif // MUELU_TENTATIVEPFACTORY_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Class that holds all level-specific information.
Factory that provides an interface for a concrete implementation of a prolongation operator.
Factory for building tentative prolongator.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void BuildPuncoupled(RCP< Matrix > A, RCP< Aggregates > aggregates, RCP< AmalgamationInfo > amalgInfo, RCP< MultiVector > fineNullspace, RCP< const Map > coarseMap, RCP< Matrix > &Ptentative, RCP< MultiVector > &coarseNullspace, const int levelID) const
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
void BuildPuncoupledBlockCrs(RCP< Matrix > A, RCP< Aggregates > aggregates, RCP< AmalgamationInfo > amalgInfo, RCP< MultiVector > fineNullspace, RCP< const Map > coarseMap, RCP< Matrix > &Ptentative, RCP< MultiVector > &coarseNullspace, const int levelID) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void BuildPcoupled(RCP< Matrix > A, RCP< Aggregates > aggregates, RCP< AmalgamationInfo > amalgInfo, RCP< MultiVector > fineNullspace, RCP< const Map > coarseMap, RCP< Matrix > &Ptentative, RCP< MultiVector > &coarseNullspace) const
Namespace for MueLu classes and methods.
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar