Loading...
Searching...
No Matches
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-2023 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 "faScalarMatrix.H"
31#include "profiling.H"
32#include "PrecisionAdaptor.H"
33
34// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35
36template<>
38(
39 const label patchI,
40 const label edgeI,
41 const direction,
42 const scalar value
43)
44{
45 const labelUList& faceLabels = psi_.mesh().boundary()[patchI].edgeFaces();
46
47 internalCoeffs_[patchI][edgeI] += diag()[faceLabels[edgeI]];
48
49 boundaryCoeffs_[patchI][edgeI] = value;
50}
51
52
53template<>
55(
56 const dictionary& solverControls
57)
58{
60 << "solving faMatrix<scalar>"
61 << endl;
62
63 const int logLevel =
64 solverControls.getOrDefault<int>
65 (
66 "log",
67 solverPerformance::debug
68 );
69
70 auto& psi = psi_.constCast();
71
72 scalarField saveDiag(diag());
73 addBoundaryDiag(diag(), 0);
74
75 scalarField totalSource(source_);
76 addBoundarySource(totalSource, 0);
77
78 // Solver call
80 (
81 psi_.name(),
82 *this,
83 boundaryCoeffs_,
84 internalCoeffs_,
85 psi_.boundaryField().scalarInterfaces(),
86 solverControls
87 )->solve(psi.primitiveFieldRef(), totalSource);
88
89 if (logLevel)
90 {
91 solverPerf.print(Info);
92 }
93
94 diag() = saveDiag;
95
96 psi.correctBoundaryConditions();
97
98 psi.mesh().data().setSolverPerformance(psi.name(), solverPerf);
99
100 return solverPerf;
101}
102
103
104template<>
106{
107 scalarField boundaryDiag(psi_.size(), Zero);
108 addBoundaryDiag(boundaryDiag, 0);
109
110 const scalarField& psif = psi_.internalField();
112 const solveScalarField& psi = tpsi();
113
115 (
117 (
118 psi,
119 source_ - boundaryDiag*psif,
120 boundaryCoeffs_,
121 psi_.boundaryField().scalarInterfaces(),
122 0
123 )
124 );
125
128
129 return tres_s;
130}
131
132
133template<>
135{
136 auto tHphi = areaScalarField::New
137 (
138 "H(" + psi_.name() + ')',
139 psi_.mesh(),
140 dimensions_/dimArea,
142 );
143 auto& Hphi = tHphi.ref();
144
145 Hphi.primitiveFieldRef() = (lduMatrix::H(psi_.primitiveField()) + source_);
146 addBoundarySource(Hphi.primitiveFieldRef());
147
148 Hphi.primitiveFieldRef() /= psi_.mesh().S();
149 Hphi.correctBoundaryConditions();
150
151 return tHphi;
152}
153
154
155// ************************************************************************* //
labelList faceLabels(nFaceLabels)
A const Field/List wrapper with possible data conversion.
static tmp< GeometricField< scalar, faPatchField, areaMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=faPatchField< scalar >::calculatedType())
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:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
SolverPerformance< Type > solve()
Solve returning the solution statistics.
const GeometricField< scalar, faPatchField, areaMesh > & psi() const
Definition faMatrix.H:348
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition faMatrix.C:619
void addBoundarySource(Field< scalar > &source, const bool couples=true) const
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
Set reference level for a component of the solution on a given patch face.
tmp< Field< Type > > residual() const
Return the matrix residual.
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
static autoPtr< solver > New(const word &solverName, const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Return a new solver of given type.
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
const scalarField & diag() const
Definition lduMatrix.C:195
tmp< Field< Type > > H(const Field< Type > &) const
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
A class for managing temporary objects.
Definition tmp.H:75
const volScalarField & psi
#define DebugInFunction
Report an information message using Foam::Info.
const dimensionSet dimArea(sqr(dimLength))
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Field< solveScalar > solveScalarField
uint8_t direction
Definition direction.H:49
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
UList< label > labelUList
A UList of labels.
Definition UList.H:75
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.