MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_BlockedCoordinatesTransferFactory_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_BLOCKEDCOORDINATESTRANSFER_FACTORY_DECL_HPP
47#define MUELU_BLOCKEDCOORDINATESTRANSFER_FACTORY_DECL_HPP
48
49#include "MueLu_ConfigDefs.hpp"
51#include "Xpetra_MultiVector_fwd.hpp"
52#include "Xpetra_MultiVectorFactory_fwd.hpp"
53#include "Xpetra_Matrix.hpp"
54
58
59namespace MueLu {
60
92 template <class Scalar = DefaultScalar,
95 class Node = DefaultNode>
97#undef MUELU_BLOCKEDCOORDINATESTRANSFERFACTORY_SHORT
99
100 public:
102
103
113
116
117 RCP<const ParameterList> GetValidParameterList() const;
118
120
122
123
129 void DeclareInput(Level &finelevel, Level &coarseLevel) const;
130
132
134
135
137 void Build(Level & fineLevel, Level &coarseLevel) const;
138
140
142
145 void AddFactory(const RCP<const FactoryBase>& factory);
146
147
149 size_t NumFactories() const { return subFactories_.size(); }
150
152 private:
154 std::vector<RCP<const FactoryBase> > subFactories_;
155
156 }; // class BlockedCoordinatesTransferFactory
157
158} // namespace MueLu
159
160#define MUELU_BLOCKEDCOORDINATESTRANSFERFACTORY_SHORT
161#endif // MUELU_BLOCKEDCOORDINATESTRANSFER_FACTORY_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Class for transferring coordinates from a finer level to a coarser one for BlockedCrsMatrices....
std::vector< RCP< const FactoryBase > > subFactories_
list of user-defined sub Factories
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
size_t NumFactories() const
Returns number of sub factories.
void DeclareInput(Level &finelevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
void AddFactory(const RCP< const FactoryBase > &factory)
Add (sub) coords factory in the end of list of factories in BlockedCoordinatesTransferFactory.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
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