Loading...
Searching...
No Matches
GAMGPreconditioner.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2013 OpenFOAM Foundation
9 Copyright (C) 2019-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "PrecisionAdaptor.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
37
38 lduMatrix::preconditioner::addsymMatrixConstructorToTable
40
41 lduMatrix::preconditioner::addasymMatrixConstructorToTable
43}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
49(
50 const lduMatrix::solver& sol,
51 const dictionary& solverControls
52)
53:
55 (
56 sol.fieldName(),
57 sol.matrix(),
58 sol.interfaceBouCoeffs(),
59 sol.interfaceIntCoeffs(),
60 sol.interfaces(),
61 solverControls
62 ),
63 lduMatrix::preconditioner(sol),
64 nVcycles_(2)
67}
68
69
70// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71
73{
75 nVcycles_ = controlDict_.getOrDefault<label>("nVcycles", 2);
76}
77
78
80(
82 const solveScalarField& rA_ss,
83 const direction cmpt
84) const
85{
86 wA = Zero;
87 solveScalarField AwA(wA.size());
88 solveScalarField finestCorrection(wA.size());
89 solveScalarField finestResidual(rA_ss);
90
91 // Create coarse grid correction fields
92 PtrList<solveScalarField> coarseCorrFields;
93
94 // Create coarse grid sources
95 PtrList<solveScalarField> coarseSources;
96
97 // Create the smoothers for all levels
99
100 // Scratch fields if processor-agglomerated coarse level meshes
101 // are bigger than original. Usually not needed
102 solveScalarField ApsiScratch;
103 solveScalarField finestCorrectionScratch;
104
105 // Initialise the above data structures
106 initVcycle
107 (
108 coarseCorrFields,
109 coarseSources,
110 smoothers,
111 ApsiScratch,
112 finestCorrectionScratch
113 );
114
115 // Adapt solveScalarField back to scalarField (as required)
117 const scalarField& rA = rA_adaptor.cref();
118
119 for (label cycle=0; cycle<nVcycles_; cycle++)
120 {
121 Vcycle
122 (
123 smoothers,
124 wA,
125 rA,
126 AwA,
127 finestCorrection,
128 finestResidual,
129
130 (ApsiScratch.size() ? ApsiScratch : AwA),
131 (
132 finestCorrectionScratch.size()
133 ? finestCorrectionScratch
134 : finestCorrection
135 ),
136
137 coarseCorrFields,
138 coarseSources,
139 cmpt
140 );
141
142 if (cycle < nVcycles_-1)
143 {
144 // Calculate finest level residual field
145 matrix_.Amul(AwA, wA, interfaceBouCoeffs_, interfaces_, cmpt);
146 finestResidual = rA_ss;
147 finestResidual -= AwA;
148 }
149 }
150}
151
152
153// ************************************************************************* //
A const Field/List wrapper with possible data conversion.
Geometric agglomerated algebraic multigrid preconditioner.
virtual void precondition(solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const
Return wA the preconditioned form of residual rA.
virtual void readControls()
Read the control parameters from the controlDict_.
GAMGPreconditioner(const lduMatrix::solver &, const dictionary &solverControls)
Construct from matrix components and preconditioner solver controls.
label nVcycles_
Number of V-cycles to perform.
Geometric agglomerated algebraic multigrid solver.
Definition GAMGSolver.H:75
friend class GAMGPreconditioner
Definition GAMGSolver.H:405
GAMGSolver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Construct from lduMatrix and solver controls.
Definition GAMGSolver.C:44
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
preconditioner(const solver &sol)
Construct for given solver.
Definition lduMatrix.H:634
Abstract base-class for lduMatrix solvers.
Definition lduMatrix.H:152
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition lduMatrix.H:338
const lduMatrix & matrix_
Definition lduMatrix.H:158
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition lduMatrix.H:343
lduInterfaceFieldPtrsList interfaces_
Definition lduMatrix.H:161
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition lduMatrix.H:159
const lduMatrix & matrix() const noexcept
Definition lduMatrix.H:328
virtual void readControls()
Read the control parameters from controlDict_.
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition lduMatrix.H:333
dictionary controlDict_
Dictionary of solution controls.
Definition lduMatrix.H:166
const word & fieldName() const noexcept
Definition lduMatrix.H:323
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition refPtrI.H:216
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
Namespace for OpenFOAM.
lduMatrix::preconditioner::addasymMatrixConstructorToTable< GAMGPreconditioner > addGAMGPreconditionerAsymMatrixConstructorToTable_
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
lduMatrix::preconditioner::addsymMatrixConstructorToTable< GAMGPreconditioner > addGAMGPreconditionerSymMatrixConstructorToTable_
Field< solveScalar > solveScalarField
uint8_t direction
Definition direction.H:49
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127