DICPreconditioner.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-2015 OpenFOAM Foundation
9 Copyright (C) 2019 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
29#include "DICPreconditioner.H"
30#include <algorithm>
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
37
38 lduMatrix::preconditioner::
39 addsymMatrixConstructorToTable<DICPreconditioner>
41}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
47(
48 const lduMatrix::solver& sol,
49 const dictionary&
50)
51:
52 lduMatrix::preconditioner(sol),
53 rD_(sol.matrix().diag().size())
54{
55 const scalarField& diag = sol.matrix().diag();
56 std::copy(diag.begin(), diag.end(), rD_.begin());
57
58 calcReciprocalD(rD_, sol.matrix());
59}
60
61
62// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63
65(
67 const lduMatrix& matrix
68)
69{
70 solveScalar* __restrict__ rDPtr = rD.begin();
71
72 const label* const __restrict__ uPtr = matrix.lduAddr().upperAddr().begin();
73 const label* const __restrict__ lPtr = matrix.lduAddr().lowerAddr().begin();
74 const scalar* const __restrict__ upperPtr = matrix.upper().begin();
75
76 // Calculate the DIC diagonal
77 const label nFaces = matrix.upper().size();
78 for (label face=0; face<nFaces; face++)
79 {
80 rDPtr[uPtr[face]] -= upperPtr[face]*upperPtr[face]/rDPtr[lPtr[face]];
81 }
82
83
84 // Calculate the reciprocal of the preconditioned diagonal
85 const label nCells = rD.size();
86
87 for (label cell=0; cell<nCells; cell++)
88 {
89 rDPtr[cell] = 1.0/rDPtr[cell];
90 }
91}
92
93
95(
97 const solveScalarField& rA,
98 const direction
99) const
100{
101 solveScalar* __restrict__ wAPtr = wA.begin();
102 const solveScalar* __restrict__ rAPtr = rA.begin();
103 const solveScalar* __restrict__ rDPtr = rD_.begin();
104
105 const label* const __restrict__ uPtr =
106 solver_.matrix().lduAddr().upperAddr().begin();
107 const label* const __restrict__ lPtr =
108 solver_.matrix().lduAddr().lowerAddr().begin();
109 const scalar* const __restrict__ upperPtr =
110 solver_.matrix().upper().begin();
111
112 const label nCells = wA.size();
113 const label nFaces = solver_.matrix().upper().size();
114 const label nFacesM1 = nFaces - 1;
115
116 for (label cell=0; cell<nCells; cell++)
117 {
118 wAPtr[cell] = rDPtr[cell]*rAPtr[cell];
119 }
120
121 for (label face=0; face<nFaces; face++)
122 {
123 wAPtr[uPtr[face]] -= rDPtr[uPtr[face]]*upperPtr[face]*wAPtr[lPtr[face]];
124 }
125
126 for (label face=nFacesM1; face>=0; face--)
127 {
128 wAPtr[lPtr[face]] -= rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[uPtr[face]];
129 }
130}
131
132
133// ************************************************************************* //
Simplified diagonal-based incomplete Cholesky preconditioner for symmetric matrices (symmetric equiva...
virtual void precondition(solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const
Return wA the preconditioned form of residual rA.
static void calcReciprocalD(solveScalarField &, const lduMatrix &)
Calculate the reciprocal of the preconditioned diagonal.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:329
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:99
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:233
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:84
scalarField & upper()
Definition: lduMatrix.C:203
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:578
scalarField & diag()
Definition: lduMatrix.C:192
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
Namespace for OpenFOAM.
uint8_t direction
Definition: direction.H:56
lduMatrix::preconditioner::addsymMatrixConstructorToTable< DICPreconditioner > addDICPreconditionerSymMatrixConstructorToTable_
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)