faScalarMatrix.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) 2016-2017 Wikki Ltd
9 Copyright (C) 2019-2021 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
27Description
28 Finite-Area scalar matrix member functions and operators
29
30\*---------------------------------------------------------------------------*/
31
32#include "faScalarMatrix.H"
34#include "PrecisionAdaptor.H"
35
36// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37
38template<>
40(
41 const label patchI,
42 const label edgeI,
43 const direction,
44 const scalar value
45)
46{
47 const labelUList& faceLabels = psi_.mesh().boundary()[patchI].edgeFaces();
48
49 internalCoeffs_[patchI][edgeI] += diag()[faceLabels[edgeI]];
50
51 boundaryCoeffs_[patchI][edgeI] = value;
52}
53
54
55template<>
57(
58 const dictionary& solverControls
59)
60{
62 << "solving faMatrix<scalar>"
63 << endl;
64
65 const int logLevel =
66 solverControls.getOrDefault<int>
67 (
68 "log",
69 solverPerformance::debug
70 );
71
72 auto& psi =
74
75 scalarField saveDiag(diag());
76 addBoundaryDiag(diag(), 0);
77
78 scalarField totalSource(source_);
79 addBoundarySource(totalSource, 0);
80
81 // Solver call
82 solverPerformance solverPerf = lduMatrix::solver::New
83 (
84 psi_.name(),
85 *this,
86 boundaryCoeffs_,
87 internalCoeffs_,
88 psi_.boundaryField().scalarInterfaces(),
89 solverControls
90 )->solve(psi.ref(), totalSource);
91
92 if (logLevel)
93 {
94 solverPerf.print(Info);
95 }
96
97 diag() = saveDiag;
98
99 psi.correctBoundaryConditions();
100
101 psi.mesh().setSolverPerformance(psi.name(), solverPerf);
102
103 return solverPerf;
104}
105
106
107template<>
109{
110 scalarField boundaryDiag(psi_.size(), Zero);
111 addBoundaryDiag(boundaryDiag, 0);
112
113 const scalarField& psif = psi_.internalField();
115 const solveScalarField& psi = tpsi();
116
118 (
119 lduMatrix::residual
120 (
121 psi,
122 source_ - boundaryDiag*psif,
123 boundaryCoeffs_,
124 psi_.boundaryField().scalarInterfaces(),
125 0
126 )
127 );
128
130 addBoundarySource(tres_s.ref());
131
132 return tres_s;
133}
134
135
136template<>
138{
140 (
142 (
144 (
145 "H("+psi_.name()+')',
146 psi_.instance(),
147 psi_.db(),
148 IOobject::NO_READ,
149 IOobject::NO_WRITE
150 ),
151 psi_.mesh(),
152 dimensions_/dimArea,
153 zeroGradientFaPatchScalarField::typeName
154 )
155 );
156 areaScalarField& Hphi = tHphi.ref();
157
158 Hphi.primitiveFieldRef() = (lduMatrix::H(psi_.primitiveField()) + source_);
159 addBoundarySource(Hphi.primitiveFieldRef());
160
161 Hphi.ref() /= psi_.mesh().S();
163
164 return tHphi;
165}
166
167
168// ************************************************************************* //
A const Field/List wrapper with possible data conversion.
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
void print(Ostream &os) const
Print summary of solver performance to the given stream.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
SolverPerformance< Type > solve()
Solve returning the solution statistics.
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition: faMatrix.C:633
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
Definition: faMatrixSolve.C:38
tmp< Field< Type > > residual() const
Return the matrix residual.
T & ref() const
Definition: refPtrI.H:203
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
const volScalarField & psi
#define DebugInFunction
Report an information message using Foam::Info.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
uint8_t direction
Definition: direction.H:56
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)