TDILUPreconditioner.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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "TDILUPreconditioner.H"
29
30// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31
32template<class Type, class DType, class LUType>
34(
36 const dictionary&
37)
38:
39 LduMatrix<Type, DType, LUType>::preconditioner(sol),
40 rD_(sol.matrix().diag())
41{
42 calcInvD(rD_, sol.matrix());
43}
44
45
46// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
47
48template<class Type, class DType, class LUType>
50(
51 Field<DType>& rD,
53)
54{
55 DType* __restrict__ rDPtr = rD.begin();
56
57 const label* const __restrict__ uPtr = matrix.lduAddr().upperAddr().begin();
58 const label* const __restrict__ lPtr = matrix.lduAddr().lowerAddr().begin();
59
60 const LUType* const __restrict__ upperPtr = matrix.upper().begin();
61 const LUType* const __restrict__ lowerPtr = matrix.lower().begin();
62
63 label nFaces = matrix.upper().size();
64 for (label face=0; face<nFaces; face++)
65 {
66 rDPtr[uPtr[face]] -=
67 dot(dot(upperPtr[face], lowerPtr[face]), inv(rDPtr[lPtr[face]]));
68 }
69
70
71 // Calculate the reciprocal of the preconditioned diagonal
72 label nCells = rD.size();
73
74 for (label cell=0; cell<nCells; cell++)
75 {
76 rDPtr[cell] = inv(rDPtr[cell]);
77 }
78}
79
80
81template<class Type, class DType, class LUType>
83(
84 Field<Type>& wA,
85 const Field<Type>& rA
86) const
87{
88 Type* __restrict__ wAPtr = wA.begin();
89 const Type* __restrict__ rAPtr = rA.begin();
90 const DType* __restrict__ rDPtr = rD_.begin();
91
92 const label* const __restrict__ uPtr =
93 this->solver_.matrix().lduAddr().upperAddr().begin();
94 const label* const __restrict__ lPtr =
95 this->solver_.matrix().lduAddr().lowerAddr().begin();
96 const label* const __restrict__ losortPtr =
97 this->solver_.matrix().lduAddr().losortAddr().begin();
98
99 const LUType* const __restrict__ upperPtr =
100 this->solver_.matrix().upper().begin();
101 const LUType* const __restrict__ lowerPtr =
102 this->solver_.matrix().lower().begin();
103
104 label nCells = wA.size();
105 label nFaces = this->solver_.matrix().upper().size();
106 label nFacesM1 = nFaces - 1;
107
108 for (label cell=0; cell<nCells; cell++)
109 {
110 wAPtr[cell] = dot(rDPtr[cell], rAPtr[cell]);
111 }
112
113
114 label sface;
115
116 for (label face=0; face<nFaces; face++)
117 {
118 sface = losortPtr[face];
119 wAPtr[uPtr[sface]] -=
120 dot(rDPtr[uPtr[sface]], dot(lowerPtr[sface], wAPtr[lPtr[sface]]));
121 }
122
123 for (label face=nFacesM1; face>=0; face--)
124 {
125 wAPtr[lPtr[face]] -=
126 dot(rDPtr[lPtr[face]], dot(upperPtr[face], wAPtr[uPtr[face]]));
127 }
128}
129
130
131template<class Type, class DType, class LUType>
133(
134 Field<Type>& wT,
135 const Field<Type>& rT
136) const
137{
138 Type* __restrict__ wTPtr = wT.begin();
139 const Type* __restrict__ rTPtr = rT.begin();
140 const DType* __restrict__ rDPtr = rD_.begin();
141
142 const label* const __restrict__ uPtr =
143 this->solver_.matrix().lduAddr().upperAddr().begin();
144 const label* const __restrict__ lPtr =
145 this->solver_.matrix().lduAddr().lowerAddr().begin();
146 const label* const __restrict__ losortPtr =
147 this->solver_.matrix().lduAddr().losortAddr().begin();
148
149 const LUType* const __restrict__ upperPtr =
150 this->solver_.matrix().upper().begin();
151 const LUType* const __restrict__ lowerPtr =
152 this->solver_.matrix().lower().begin();
153
154 label nCells = wT.size();
155 label nFaces = this->solver_.matrix().upper().size();
156 label nFacesM1 = nFaces - 1;
157
158 for (label cell=0; cell<nCells; cell++)
159 {
160 wTPtr[cell] = dot(rDPtr[cell], rTPtr[cell]);
161 }
162
163 for (label face=0; face<nFaces; face++)
164 {
165 wTPtr[uPtr[face]] -=
166 dot(rDPtr[uPtr[face]], dot(upperPtr[face], wTPtr[lPtr[face]]));
167 }
168
169
170 label sface;
171
172 for (label face=nFacesM1; face>=0; face--)
173 {
174 sface = losortPtr[face];
175 wTPtr[lPtr[sface]] -=
176 dot(rDPtr[lPtr[sface]], dot(lowerPtr[sface], wTPtr[uPtr[sface]]));
177 }
178}
179
180
181// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
Abstract base-class for LduMatrix solvers.
Definition: LduMatrix.H:116
const LduMatrix< Type, DType, LUType > & matrix() const noexcept
Definition: LduMatrix.H:239
LduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: LduMatrix.H:88
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: LduMatrix.H:498
Field< LUType > & upper()
Definition: LduMatrix.C:204
Field< LUType > & lower()
Definition: LduMatrix.C:227
Simplified diagonal-based incomplete LU preconditioner for asymmetric matrices.
static void calcInvD(Field< DType > &rD, const LduMatrix< Type, DType, LUType > &matrix)
Calculate the reciprocal of the preconditioned diagonal.
virtual void precondition(Field< Type > &wA, const Field< Type > &rA) const
Return wA the preconditioned form of residual rA.
virtual void preconditionT(Field< Type > &wT, const Field< Type > &rT) const
Return wT the transpose-matrix preconditioned form of.
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.
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)