MueLu
Version of the Day
Toggle main menu visibility
Loading...
Searching...
No Matches
MueLu_SemiCoarsenPFactory_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_SEMICOARSENPFACTORY_DECL_HPP
47
#define MUELU_SEMICOARSENPFACTORY_DECL_HPP
48
49
#include <Xpetra_Matrix_fwd.hpp>
50
51
#include "
MueLu_ConfigDefs.hpp
"
52
#include "MueLu_PFactory.hpp"
53
#include "
MueLu_SemiCoarsenPFactory_fwd.hpp
"
54
55
#include "
MueLu_Level_fwd.hpp
"
56
57
#define VERTICAL 1
58
#define HORIZONTAL 2
59
#define GRID_SUPPLIED -1
60
#define NUM_ZPTS 0
61
#define ORIENTATION 1
62
63
namespace
MueLu
{
64
65
107
template
<
class
Scalar
=
DefaultScalar
,
108
class
LocalOrdinal
=
DefaultLocalOrdinal
,
109
class
GlobalOrdinal
=
DefaultGlobalOrdinal
,
110
class
Node
=
DefaultNode
>
111
class
SemiCoarsenPFactory
:
public
PFactory
{
112
#undef MUELU_SEMICOARSENPFACTORY_SHORT
113
#include "
MueLu_UseShortNames.hpp
"
114
115
public
:
117
118
120
SemiCoarsenPFactory
() :
bTransferCoordinates_
(false) { }
121
123
virtual
~SemiCoarsenPFactory
() { }
125
126
RCP<const ParameterList>
GetValidParameterList
()
const
;
127
129
130
131
void
DeclareInput
(
Level
& fineLevel,
Level
& coarseLevel)
const
;
132
134
136
137
138
void
Build
(
Level
& fineLevel,
Level
& coarseLevel)
const
;
139
void
BuildP
(
Level
& fineLevel,
Level
& coarseLevel)
const
;
140
142
143
private
:
144
LO
FindCpts
(LO
const
PtsPerLine, LO
const
CoarsenRate, LO
const
Thin, LO **LayerCpts)
const
;
145
LO
MakeSemiCoarsenP
(LO
const
Ntotal, LO
const
nz, LO
const
CoarsenRate, LO
const
LayerId[],
146
LO
const
VertLineId[], LO
const
DofsPerNode, RCP<Matrix>& Amat,
147
RCP<Matrix>& P, RCP<const Map>& coarseMap,
148
const
RCP<MultiVector> fineNullspace, RCP<MultiVector>& coarseNullspace, RCP<Matrix>& R,
bool
buildRestriction,
bool
doLinear)
const
;
149
void
RevertToPieceWiseConstant
( RCP<Matrix>& P, LO BlkSize)
const
;
150
151
152
mutable
bool
bTransferCoordinates_
;
//< boolean which is true if coordinate information is available to be transferred to coarse coordinate information
153
};
//class SemiCoarsenPFactory
154
155
}
//namespace MueLu
156
157
#define MUELU_SEMICOARSENPFACTORY_SHORT
158
#endif
// MUELU_SEMICOARSENPFACTORY_DECL_HPP
MueLu_ConfigDefs.hpp
MueLu_Level_fwd.hpp
MueLu_SemiCoarsenPFactory_fwd.hpp
LocalOrdinal
MueLu::DefaultLocalOrdinal LocalOrdinal
Definition
MueLu_UseDefaultTypes.hpp:50
Scalar
MueLu::DefaultScalar Scalar
Definition
MueLu_UseDefaultTypes.hpp:49
GlobalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Definition
MueLu_UseDefaultTypes.hpp:51
Node
MueLu::DefaultNode Node
Definition
MueLu_UseDefaultTypes.hpp:52
MueLu_UseShortNames.hpp
MueLu::Level
Class that holds all level-specific information.
Definition
MueLu_Level.hpp:99
MueLu::PFactory::PFactory
PFactory()
Constructor.
Definition
MueLu_PFactory_decl.hpp:75
MueLu::SemiCoarsenPFactory::FindCpts
LO FindCpts(LO const PtsPerLine, LO const CoarsenRate, LO const Thin, LO **LayerCpts) const
Definition
MueLu_SemiCoarsenPFactory_def.hpp:280
MueLu::SemiCoarsenPFactory::SemiCoarsenPFactory
SemiCoarsenPFactory()
Constructor.
Definition
MueLu_SemiCoarsenPFactory_decl.hpp:120
MueLu::SemiCoarsenPFactory::DeclareInput
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
Definition
MueLu_SemiCoarsenPFactory_def.hpp:89
MueLu::SemiCoarsenPFactory::GetValidParameterList
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Definition
MueLu_SemiCoarsenPFactory_def.hpp:68
MueLu::SemiCoarsenPFactory::MakeSemiCoarsenP
LO MakeSemiCoarsenP(LO const Ntotal, LO const nz, LO const CoarsenRate, LO const LayerId[], LO const VertLineId[], LO const DofsPerNode, RCP< Matrix > &Amat, RCP< Matrix > &P, RCP< const Map > &coarseMap, const RCP< MultiVector > fineNullspace, RCP< MultiVector > &coarseNullspace, RCP< Matrix > &R, bool buildRestriction, bool doLinear) const
Definition
MueLu_SemiCoarsenPFactory_def.hpp:336
MueLu::SemiCoarsenPFactory::BuildP
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
Definition
MueLu_SemiCoarsenPFactory_def.hpp:123
MueLu::SemiCoarsenPFactory::~SemiCoarsenPFactory
virtual ~SemiCoarsenPFactory()
Destructor.
Definition
MueLu_SemiCoarsenPFactory_decl.hpp:123
MueLu::SemiCoarsenPFactory::Build
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
Definition
MueLu_SemiCoarsenPFactory_def.hpp:118
MueLu::SemiCoarsenPFactory::RevertToPieceWiseConstant
void RevertToPieceWiseConstant(RCP< Matrix > &P, LO BlkSize) const
Definition
MueLu_SemiCoarsenPFactory_def.hpp:960
MueLu::SemiCoarsenPFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node >::bTransferCoordinates_
bool bTransferCoordinates_
Definition
MueLu_SemiCoarsenPFactory_decl.hpp:152
MueLu
Namespace for MueLu classes and methods.
Definition
MueLu_BrickAggregationFactory_decl.hpp:78
MueLu::DefaultNode
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Definition
MueLu_Details_DefaultTypes.hpp:70
MueLu::DefaultGlobalOrdinal
int DefaultGlobalOrdinal
Definition
MueLu_Details_DefaultTypes.hpp:67
MueLu::DefaultLocalOrdinal
int DefaultLocalOrdinal
Definition
MueLu_Details_DefaultTypes.hpp:60
MueLu::DefaultScalar
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar
Definition
MueLu_Details_DefaultTypes.hpp:58
src
Transfers
SemiCoarsen
MueLu_SemiCoarsenPFactory_decl.hpp
Generated by
1.17.0