Nektar++
Preconditioner.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Preconditioner.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Preconditioner header
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_LIB_MULTIREGIONS_PRECONDITIONER_H
36#define NEKTAR_LIB_MULTIREGIONS_PRECONDITIONER_H
37
44
45#include <memory>
46
47namespace Nektar
48{
49namespace MultiRegions
50{
51class AssemblyMap;
52typedef std::shared_ptr<AssemblyMap> AssemblyMapSharedPtr;
53
54class Preconditioner;
55typedef std::shared_ptr<Preconditioner> PreconditionerSharedPtr;
56
58
59typedef LibUtilities::NekFactory<std::string, Preconditioner,
60 const std::shared_ptr<GlobalLinSys> &,
61 const std::shared_ptr<AssemblyMap> &>
64
66{
67public:
69 const std::shared_ptr<GlobalLinSys> &plinsys,
70 const AssemblyMapSharedPtr &pLocToGloMap);
71
74 {
75 }
76
77 inline void DoPreconditioner(const Array<OneD, NekDouble> &pInput,
79 const bool &IsLocal = false);
80
81 void DoAssembleLoc(const Array<OneD, NekDouble> &pInput,
82 Array<OneD, NekDouble> &pOutput, const bool &ZeroDir);
83
85 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
86 const Array<OneD, NekDouble> &pNonVertOutput,
88
90
92
94 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput);
95
97 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput);
98
99 inline void BuildPreconditioner();
100
101 inline void InitObject();
102
104
106
108
110
112
114
116 const;
117
119 int offset, int bndoffset, const std::shared_ptr<DNekScalMat> &loc_mat);
120
121protected:
122 const std::weak_ptr<GlobalLinSys> m_linsys;
123 std::string m_preconType;
125 std::weak_ptr<AssemblyMap> m_locToGloMap;
127
129 int offset, int bndoffset, const std::shared_ptr<DNekScalMat> &loc_mat);
130
131 virtual void v_InitObject();
132
133 virtual void v_DoPreconditioner(const Array<OneD, NekDouble> &pInput,
134 Array<OneD, NekDouble> &pOutput,
135 const bool &isLocal = false);
136
138 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
139 const Array<OneD, NekDouble> &pNonVertOutput,
140 Array<OneD, NekDouble> &pVertForce);
141
143
145 Array<OneD, NekDouble> &pInOut);
146
148 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput);
149
151 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput);
152
153 virtual void v_BuildPreconditioner();
154
155private:
157
158 static std::string lookupIds[];
159 static std::string def;
160};
161typedef std::shared_ptr<Preconditioner> PreconditionerSharedPtr;
162
163/**
164 *
165 */
167{
168 v_InitObject();
169}
170
171/**
172 *
173 */
175 int offset, int bndoffset, const std::shared_ptr<DNekScalMat> &loc_mat)
176{
177 return v_TransformedSchurCompl(offset, bndoffset, loc_mat);
178}
179
180/**
181 *
182 */
184 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
185 const bool &ZeroDir)
186{
187 v_DoPreconditioner(pInput, pOutput, ZeroDir);
188}
189
190/**
191 *
192 */
194 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
195 const Array<OneD, NekDouble> &pNonVertOutput,
196 Array<OneD, NekDouble> &pVertForce)
197{
198 v_DoPreconditionerWithNonVertOutput(pInput, pOutput, pNonVertOutput,
199 pVertForce);
200}
201
202/**
203 *
204 */
207{
209}
210
211/**
212 *
213 */
216{
218}
219
220/**
221 *
222 */
224 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
225{
226 v_DoTransformCoeffsToLowEnergy(pInput, pOutput);
227}
228
229/**
230 *
231 */
233 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
234{
235 v_DoTransformBasisFromLowEnergy(pInput, pOutput);
236}
237
238/**
239 *
240 */
242{
244}
245} // namespace MultiRegions
246} // namespace Nektar
247
248#endif
#define MULTI_REGIONS_EXPORT
Provides a generic Factory class.
Definition: NekFactory.hpp:105
const DNekScalBlkMatSharedPtr & GetBlockTransformationMatrix() const
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false)
Apply a preconditioner to the conjugate gradient method.
void DoTransformCoeffsToLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
virtual void v_DoTransformBasisToLowEnergy(Array< OneD, NekDouble > &pInOut)
Transform from original basis to low energy basis.
LibUtilities::CommSharedPtr m_comm
const std::weak_ptr< GlobalLinSys > m_linsys
const DNekScalBlkMatSharedPtr & GetBlockTransformedSchurCompl() const
virtual void v_DoPreconditionerWithNonVertOutput(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce)
Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of free...
const DNekScalBlkMatSharedPtr & GetBlockTransposedTransformationMatrix() const
Array< OneD, NekDouble > AssembleStaticCondGlobalDiagonals()
Performs global assembly of diagonal entries to global Schur complement matrix.
void DoAssembleLoc(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &ZeroDir)
Apply an assembly and scatter back to lcoal array.
std::weak_ptr< AssemblyMap > m_locToGloMap
virtual void v_DoTransformBasisFromLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block transposed inverse transformation matrix.
void DoTransformBasisToLowEnergy(Array< OneD, NekDouble > &pInOut)
const DNekScalBlkMatSharedPtr & GetBlockCMatrix() const
void DoPreconditionerWithNonVertOutput(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce=NullNekDouble1DArray)
DNekScalMatSharedPtr TransformedSchurCompl(int offset, int bndoffset, const std::shared_ptr< DNekScalMat > &loc_mat)
Preconditioner(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
void DoTransformBasisFromLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
virtual void v_DoTransformCoeffsFromLowEnergy(Array< OneD, NekDouble > &pInOut)
Transform from low energy coeffs to orignal basis.
virtual DNekScalMatSharedPtr v_TransformedSchurCompl(int offset, int bndoffset, const std::shared_ptr< DNekScalMat > &loc_mat)
Get block elemental transposed transformation matrix .
void DoTransformCoeffsFromLowEnergy(Array< OneD, NekDouble > &pInOut)
const DNekScalBlkMatSharedPtr & GetBlockSchurCompl() const
const DNekScalBlkMatSharedPtr & GetBlockInvDMatrix() const
void DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &IsLocal=false)
virtual void v_DoTransformCoeffsToLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block inverse transformation matrix.
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
static PreconditionerSharedPtr NullPreconditionerSharedPtr
LibUtilities::NekFactory< std::string, Preconditioner, const std::shared_ptr< GlobalLinSys > &, const std::shared_ptr< AssemblyMap > & > PreconFactory
std::shared_ptr< Preconditioner > PreconditionerSharedPtr
Definition: GlobalLinSys.h:60
PreconFactory & GetPreconFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75