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 -------------------------------------------------------------------------------
11 License
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 
36 template<>
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 
53 template<>
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",
68  );
69 
70  auto& psi =
72 
73  scalarField saveDiag(diag());
74  addBoundaryDiag(diag(), 0);
75 
76  scalarField totalSource(source_);
77  addBoundarySource(totalSource, 0);
78 
79  // Solver call
81  (
82  psi_.name(),
83  *this,
84  boundaryCoeffs_,
85  internalCoeffs_,
86  psi_.boundaryField().scalarInterfaces(),
87  solverControls
88  )->solve(psi.ref(), totalSource);
89 
90  if (logLevel)
91  {
92  solverPerf.print(Info);
93  }
94 
95  diag() = saveDiag;
96 
97  psi.correctBoundaryConditions();
98 
99  psi.mesh().data().setSolverPerformance(psi.name(), solverPerf);
100 
101  return solverPerf;
102 }
103 
104 
105 template<>
107 {
108  scalarField boundaryDiag(psi_.size(), Zero);
109  addBoundaryDiag(boundaryDiag, 0);
110 
111  const scalarField& psif = psi_.internalField();
113  const solveScalarField& psi = tpsi();
114 
116  (
117  lduMatrix::residual
118  (
119  psi,
120  source_ - boundaryDiag*psif,
121  boundaryCoeffs_,
122  psi_.boundaryField().scalarInterfaces(),
123  0
124  )
125  );
126 
128  addBoundarySource(tres_s.ref());
129 
130  return tres_s;
131 }
132 
133 
134 template<>
136 {
137  auto tHphi = areaScalarField::New
138  (
139  "H(" + psi_.name() + ')',
140  psi_.mesh(),
141  dimensions_/dimArea,
143  );
144  auto& Hphi = tHphi.ref();
145 
146  Hphi.primitiveFieldRef() = (lduMatrix::H(psi_.primitiveField()) + source_);
147  addBoundarySource(Hphi.primitiveFieldRef());
148 
149  Hphi.ref() /= psi_.mesh().S();
150  Hphi.correctBoundaryConditions();
151 
152  return tHphi;
153 }
154 
155 
156 // ************************************************************************* //
uint8_t direction
Definition: direction.H:46
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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.
Definition: faMatrixSolve.C:32
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
labelList faceLabels(nFaceLabels)
CEqn solve()
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
SolverPerformance< Type > solve()
Solve returning the solution statistics.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
#define DebugInFunction
Report an information message using Foam::Info.
A const Field/List wrapper with possible data conversion.
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition: faMatrix.C:634
int debug
Static debugging option.
Container< Type > & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: refPtrI.H:230
void print(Ostream &os) const
Print summary of solver performance to the given stream.
tmp< Field< Type > > residual() const
Return the matrix residual.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const volScalarField & psi
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...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const word extrapolatedCalculatedType
A combined zero-gradient and calculated patch field type.
volScalarField H(IOobject("H", runTime.timeName(), mesh.thisDb(), IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127