Files
openfoam/src/OpenFOAM/matrices/lduMatrix/preconditioners/DICPreconditioner/DICPreconditioner.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

184 lines
5.8 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 "DICPreconditioner.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(DICPreconditioner, 0);
lduMatrix::preconditioner::
addsymMatrixConstructorToTable<DICPreconditioner>
addDICPreconditionerSymMatrixConstructorToTable_;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::DICPreconditioner::DICPreconditioner
(
const lduMatrix::solver& sol,
const dictionary&
)
:
lduMatrix::preconditioner(sol),
rD_(sol.matrix().diag())
{
calcReciprocalD(rD_, sol.matrix());
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::DICPreconditioner::calcReciprocalD
(
scalarField& rD,
const lduMatrix& matrix
)
{
scalar* __restrict__ rDPtr = rD.begin();
const label* const __restrict__ uPtr = matrix.lduAddr().upperAddr().begin();
const label* const __restrict__ lPtr = matrix.lduAddr().lowerAddr().begin();
const scalar* const __restrict__ upperPtr = matrix.upper().begin();
// Calculate the DIC diagonal
register const label nFaces = matrix.upper().size();
for (register label face=0; face<nFaces; face++)
{
#ifdef ICC_IA64_PREFETCH
__builtin_prefetch (&uPtr[face+96],0,0);
__builtin_prefetch (&lPtr[face+96],0,0);
__builtin_prefetch (&upperPtr[face+96],0,1);
__builtin_prefetch (&rDPtr[lPtr[face+24]],0,1);
__builtin_prefetch (&rDPtr[uPtr[face+24]],1,1);
#endif
rDPtr[uPtr[face]] -= upperPtr[face]*upperPtr[face]/rDPtr[lPtr[face]];
}
// Calculate the reciprocal of the preconditioned diagonal
register const label nCells = rD.size();
#ifdef ICC_IA64_PREFETCH
#pragma ivdep
#endif
for (register label cell=0; cell<nCells; cell++)
{
#ifdef ICC_IA64_PREFETCH
__builtin_prefetch (&rDPtr[cell+96],0,1);
#endif
rDPtr[cell] = 1.0/rDPtr[cell];
}
}
void Foam::DICPreconditioner::precondition
(
scalarField& wA,
const scalarField& rA,
const direction
) const
{
scalar* __restrict__ wAPtr = wA.begin();
const scalar* __restrict__ rAPtr = rA.begin();
const scalar* __restrict__ rDPtr = rD_.begin();
const label* const __restrict__ uPtr =
solver_.matrix().lduAddr().upperAddr().begin();
const label* const __restrict__ lPtr =
solver_.matrix().lduAddr().lowerAddr().begin();
const scalar* const __restrict__ upperPtr = solver_.matrix().upper().begin();
register label nCells = wA.size();
register label nFaces = solver_.matrix().upper().size();
register label nFacesM1 = nFaces - 1;
#ifdef ICC_IA64_PREFETCH
#pragma ivdep
#endif
for (register label cell=0; cell<nCells; cell++)
{
#ifdef ICC_IA64_PREFETCH
__builtin_prefetch (&wAPtr[cell+96],0,1);
__builtin_prefetch (&rDPtr[cell+96],0,1);
__builtin_prefetch (&rAPtr[cell+96],0,1);
#endif
wAPtr[cell] = rDPtr[cell]*rAPtr[cell];
}
#ifdef ICC_IA64_PREFETCH
#pragma noprefetch uPtr,lPtr,upperPtr,rDPtr,wAPtr
#pragma nounroll
#endif
for (register label face=0; face<nFaces; face++)
{
#ifdef ICC_IA64_PREFETCH
__builtin_prefetch (&uPtr[face+96],0,0);
__builtin_prefetch (&lPtr[face+96],0,0);
__builtin_prefetch (&upperPtr[face+96],0,0);
__builtin_prefetch (&rDPtr[uPtr[face+32]],0,1);
__builtin_prefetch (&wAPtr[uPtr[face+32]],0,1);
__builtin_prefetch (&wAPtr[lPtr[face+32]],0,1);
#endif
wAPtr[uPtr[face]] -= rDPtr[uPtr[face]]*upperPtr[face]*wAPtr[lPtr[face]];
}
#ifdef ICC_IA64_PREFETCH
#pragma noprefetch uPtr,lPtr,rDPtr,wAPtr
#pragma nounroll
#endif
for (register label face=nFacesM1; face>=0; face--)
{
#ifdef ICC_IA64_PREFETCH
__builtin_prefetch (&uPtr[face-95],0,0);
__builtin_prefetch (&lPtr[face-95],0,0);
__builtin_prefetch (&rDPtr[lPtr[face-16]],0,1);
__builtin_prefetch (&wAPtr[lPtr[face-16]],0,1);
__builtin_prefetch (&wAPtr[uPtr[face-16]],0,1);
__builtin_prefetch (&rDPtr[lPtr[face-24]],0,1);
__builtin_prefetch (&wAPtr[lPtr[face-24]],0,1);
__builtin_prefetch (&wAPtr[uPtr[face-24]],0,1);
__builtin_prefetch (&rDPtr[lPtr[face-32]],0,1);
__builtin_prefetch (&wAPtr[lPtr[face-32]],0,1);
__builtin_prefetch (&wAPtr[uPtr[face-32]],0,1);
#endif
wAPtr[lPtr[face]] -= rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[uPtr[face]];
}
}
// ************************************************************************* //