Files
openfoam/src/OpenFOAM/matrices/lduMatrix/preconditioners/GAMGPreconditioner/GAMGPreconditioner.C
Mark Olesen 4f1bb8345f lduMatrix now takes a dictionary instead of an Istream for the solver controls
- can now use dictionary substitutions and regular expressions in
    system/fvSolution

  - foamUpgradeFvSolution application to convert system/fvSolution
    (with -test option)

motion solver syntax left as-is.
2008-12-08 17:22:01 +01:00

133 lines
3.7 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "GAMGPreconditioner.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(GAMGPreconditioner, 0);
lduMatrix::preconditioner::addsymMatrixConstructorToTable
<GAMGPreconditioner> addGAMGPreconditionerSymMatrixConstructorToTable_;
lduMatrix::preconditioner::addasymMatrixConstructorToTable
<GAMGPreconditioner> addGAMGPreconditionerAsymMatrixConstructorToTable_;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::GAMGPreconditioner::GAMGPreconditioner
(
const lduMatrix::solver& sol,
const dictionary& solverControls
)
:
GAMGSolver
(
sol.fieldName(),
sol.matrix(),
sol.interfaceBouCoeffs(),
sol.interfaceIntCoeffs(),
sol.interfaces(),
solverControls
),
lduMatrix::preconditioner(sol),
nVcycles_(2)
{
readControls();
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::GAMGPreconditioner::~GAMGPreconditioner()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::GAMGPreconditioner::readControls()
{
GAMGSolver::readControls();
nVcycles_ = controlDict_.lookupOrDefault<label>("nVcycles", 2);
}
void Foam::GAMGPreconditioner::precondition
(
scalarField& wA,
const scalarField& rA,
const direction cmpt
) const
{
wA = 0.0;
scalarField AwA(wA.size());
scalarField finestCorrection(wA.size());
scalarField finestResidual(rA);
// Create coarse grid correction fields
PtrList<scalarField> coarseCorrFields;
// Create coarse grid sources
PtrList<scalarField> coarseSources;
// Create the smoothers for all levels
PtrList<lduMatrix::smoother> smoothers;
// Initialise the above data structures
initVcycle(coarseCorrFields, coarseSources, smoothers);
for (label cycle=0; cycle<nVcycles_; cycle++)
{
Vcycle
(
smoothers,
wA,
rA,
AwA,
finestCorrection,
finestResidual,
coarseCorrFields,
coarseSources,
cmpt
);
if (cycle < nVcycles_-1)
{
// Calculate finest level residual field
matrix_.Amul(AwA, wA, interfaceBouCoeffs_, interfaces_, cmpt);
finestResidual = rA;
finestResidual -= AwA;
}
}
}
// ************************************************************************* //