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 -------------------------------------------------------------------------------
11 License
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 
29 #include "GAMGPreconditioner.H"
30 #include "PrecisionAdaptor.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
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)
65 {
66  readControls();
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
71 
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 {
80  GAMGSolver::readControls();
81  nVcycles_ = controlDict_.getOrDefault<label>("nVcycles", 2);
82 }
83 
84 
86 (
87  solveScalarField& wA,
88  const solveScalarField& rA_ss,
89  const direction cmpt
90 ) const
91 {
92  wA = 0.0;
93  solveScalarField AwA(wA.size());
94  solveScalarField finestCorrection(wA.size());
95  solveScalarField finestResidual(rA_ss);
96 
97  // Create coarse grid correction fields
98  PtrList<solveScalarField> coarseCorrFields;
99 
100  // Create coarse grid sources
101  PtrList<solveScalarField> coarseSources;
102 
103  // Create the smoothers for all levels
105 
106  // Scratch fields if processor-agglomerated coarse level meshes
107  // are bigger than original. Usually not needed
108  solveScalarField ApsiScratch;
109  solveScalarField finestCorrectionScratch;
110 
111  // Initialise the above data structures
112  initVcycle
113  (
114  coarseCorrFields,
115  coarseSources,
116  smoothers,
117  ApsiScratch,
118  finestCorrectionScratch
119  );
120 
121  // Adapt solveScalarField back to scalarField (as required)
123  const scalarField& rA = rA_adaptor.cref();
124 
125  for (label cycle=0; cycle<nVcycles_; cycle++)
126  {
127  Vcycle
128  (
129  smoothers,
130  wA,
131  rA,
132  AwA,
133  finestCorrection,
134  finestResidual,
135 
136  (ApsiScratch.size() ? ApsiScratch : AwA),
137  (
138  finestCorrectionScratch.size()
139  ? finestCorrectionScratch
140  : finestCorrection
141  ),
142 
143  coarseCorrFields,
144  coarseSources,
145  cmpt
146  );
147 
148  if (cycle < nVcycles_-1)
149  {
150  // Calculate finest level residual field
151  matrix_.Amul(AwA, wA, interfaceBouCoeffs_, interfaces_, cmpt);
152  finestResidual = rA_ss;
153  finestResidual -= AwA;
154  }
155  }
156 }
157 
158 
159 // ************************************************************************* //
Foam::lduMatrix::solver::interfaces
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:244
Foam::GAMGSolver::GAMGPreconditioner
friend class GAMGPreconditioner
Definition: GAMGSolver.H:325
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:79
Foam::addGAMGPreconditionerSymMatrixConstructorToTable_
lduMatrix::preconditioner::addsymMatrixConstructorToTable< GAMGPreconditioner > addGAMGPreconditionerSymMatrixConstructorToTable_
Definition: GAMGPreconditioner.C:32
PrecisionAdaptor.H
Foam::lduMatrix::solver
Definition: lduMatrix.H:94
Foam::lduMatrix::solver::fieldName
const word & fieldName() const noexcept
Definition: lduMatrix.H:224
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::refPtr< Field< Type > >::cref
const T & cref() const
Definition: refPtrI.H:182
Foam::lduMatrix::solver::interfaceBouCoeffs
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition: lduMatrix.H:234
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
Foam::GAMGPreconditioner::readControls
virtual void readControls()
Definition: GAMGPreconditioner.C:71
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
Foam::GAMGPreconditioner::~GAMGPreconditioner
virtual ~GAMGPreconditioner()
Definition: GAMGPreconditioner.C:65
Foam::GAMGPreconditioner
Geometric agglomerated algebraic multigrid preconditioner.
Definition: GAMGPreconditioner.H:53
Foam
Definition: atmBoundaryLayer.C:26
Foam::GAMGPreconditioner::precondition
virtual void precondition(solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const
Definition: GAMGPreconditioner.C:79
Foam::GAMGSolver
Geometric agglomerated algebraic multigrid solver.
Definition: GAMGSolver.H:72
Foam::lduMatrix::solver::interfaceIntCoeffs
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition: lduMatrix.H:239
Foam::ConstPrecisionAdaptor
Definition: PrecisionAdaptor.H:53
Foam::direction
uint8_t direction
Definition: direction.H:46
Foam::lduMatrix::preconditioner
Definition: lduMatrix.H:429
Foam::addGAMGPreconditionerAsymMatrixConstructorToTable_
lduMatrix::preconditioner::addasymMatrixConstructorToTable< GAMGPreconditioner > addGAMGPreconditionerAsymMatrixConstructorToTable_
Definition: GAMGPreconditioner.C:35
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
GAMGPreconditioner.H
Foam::lduMatrix::solver::matrix
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:229