Update iluSmoother for symmetric and asymmetric matrices

This commit is contained in:
Hrvoje Jasak 2017-06-19 17:18:50 +01:00
parent 0cbbf0d443
commit 1354e6971c
2 changed files with 41 additions and 16 deletions

View file

@ -25,7 +25,7 @@ Class
iluSmoother
Description
Symmetric Gauss-Seidel smoother
ILU smoother
Author
Hrvoje Jasak, Wikki Ltd. All rights reserved.
@ -33,6 +33,8 @@ Author
\*---------------------------------------------------------------------------*/
#include "iluSmoother.H"
#include "CholeskyPrecon.H"
#include "ILU0.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
@ -66,16 +68,37 @@ Foam::iluSmoother::iluSmoother
coupleIntCoeffs,
interfaces
),
precon_
(
matrix,
coupleBouCoeffs,
coupleIntCoeffs,
interfaces
),
preconPtr_(),
xCorr_(matrix.lduAddr().size()),
residual_(matrix.lduAddr().size())
{}
{
if (matrix.symmetric())
{
preconPtr_.set
(
new CholeskyPrecon
(
matrix,
coupleBouCoeffs,
coupleIntCoeffs,
interfaces
)
);
}
else if (matrix.asymmetric())
{
preconPtr_.set
(
new ILU0
(
matrix,
coupleBouCoeffs,
coupleIntCoeffs,
interfaces
)
);
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
@ -106,7 +129,10 @@ void Foam::iluSmoother::smooth
residual_[i] = b[i] - residual_[i];
}
precon_.precondition(xCorr_, residual_, cmpt);
if (!matrix_.diagonal())
{
preconPtr_->precondition(xCorr_, residual_, cmpt);
}
// Add correction to x
x += xCorr_;

View file

@ -38,7 +38,7 @@ SourceFiles
#ifndef iluSmoother_H
#define iluSmoother_H
#include "CholeskyPrecon.H"
#include "lduMatrix.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -56,7 +56,7 @@ class iluSmoother
// Private Data
//- Cholesky preconditioner
CholeskyPrecon precon_;
autoPtr<lduPreconditioner> preconPtr_;
//- Correction array
mutable scalarField xCorr_;
@ -92,10 +92,9 @@ public:
);
// Destructor
virtual ~iluSmoother()
{}
//- Destructor
virtual ~iluSmoother()
{}
// Member Functions