fvScalarMatrix.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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2019 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 "fvScalarMatrix.H"
31 #include "profiling.H"
32 #include "PrecisionAdaptor.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<>
38 (
39  const label patchi,
40  const label facei,
41  const direction,
42  const scalar value
43 )
44 {
45  if (psi_.needReference())
46  {
47  if (Pstream::master())
48  {
49  internalCoeffs_[patchi][facei] +=
50  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
51 
52  boundaryCoeffs_[patchi][facei] +=
53  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]
54  *value;
55  }
56  }
57 }
58 
59 
60 template<>
63 (
64  const dictionary& solverControls
65 )
66 {
68  if (psi_.mesh().name() != polyMesh::defaultRegion)
69  {
70  regionName = psi_.mesh().name() + "::";
71  }
72  addProfiling(solve, "fvMatrix::solve." + regionName + psi_.name());
73 
74  if (debug)
75  {
76  Info.masterStream(this->mesh().comm())
77  << "fvMatrix<scalar>::solver(const dictionary& solverControls) : "
78  "solver for fvMatrix<scalar>"
79  << endl;
80  }
81 
82  scalarField saveDiag(diag());
83  addBoundaryDiag(diag(), 0);
84 
86  (
88  (
89  *this,
91  (
92  psi_.name(),
93  *this,
94  boundaryCoeffs_,
95  internalCoeffs_,
96  psi_.boundaryField().scalarInterfaces(),
97  solverControls
98  )
99  )
100  );
101 
102  diag() = saveDiag;
103 
104  return solverPtr;
105 }
106 
107 
108 template<>
110 (
111  const dictionary& solverControls
112 )
113 {
116  (fvMat_.psi());
117 
118  scalarField saveDiag(fvMat_.diag());
119  fvMat_.addBoundaryDiag(fvMat_.diag(), 0);
120 
121  scalarField totalSource(fvMat_.source());
122  fvMat_.addBoundarySource(totalSource, false);
123 
124  // Assign new solver controls
125  solver_->read(solverControls);
126 
127  solverPerformance solverPerf = solver_->solve
128  (
129  psi.primitiveFieldRef(),
130  totalSource
131  );
132 
134  {
135  solverPerf.print(Info.masterStream(fvMat_.mesh().comm()));
136  }
137 
138  fvMat_.diag() = saveDiag;
139 
140  psi.correctBoundaryConditions();
141 
142  psi.mesh().setSolverPerformance(psi.name(), solverPerf);
143 
144  return solverPerf;
145 }
146 
147 
148 template<>
150 (
151  const dictionary& solverControls
152 )
153 {
154  if (debug)
155  {
156  Info.masterStream(this->mesh().comm())
157  << "fvMatrix<scalar>::solveSegregated"
158  "(const dictionary& solverControls) : "
159  "solving fvMatrix<scalar>"
160  << endl;
161  }
162 
165 
166  scalarField saveDiag(diag());
167  addBoundaryDiag(diag(), 0);
168 
169  scalarField totalSource(source_);
170  addBoundarySource(totalSource, false);
171 
172  // Solver call
174  (
175  psi.name(),
176  *this,
177  boundaryCoeffs_,
178  internalCoeffs_,
179  psi_.boundaryField().scalarInterfaces(),
180  solverControls
181  )->solve(psi.primitiveFieldRef(), totalSource);
182 
184  {
185  solverPerf.print(Info.masterStream(mesh().comm()));
186  }
187 
188  diag() = saveDiag;
189 
190  psi.correctBoundaryConditions();
191 
192  psi.mesh().setSolverPerformance(psi.name(), solverPerf);
193 
194  return solverPerf;
195 }
196 
197 
198 template<>
200 {
201  scalarField boundaryDiag(psi_.size(), Zero);
202  addBoundaryDiag(boundaryDiag, 0);
203 
204  const scalarField& psif = psi_.primitiveField();
206  const solveScalarField& psi = tpsi();
207 
209  (
210  lduMatrix::residual
211  (
212  psi,
213  source_ - boundaryDiag*psif,
214  boundaryCoeffs_,
215  psi_.boundaryField().scalarInterfaces(),
216  0
217  )
218  );
219 
221  addBoundarySource(tres_s.ref());
222 
223  return tres_s;
224 }
225 
226 
227 template<>
229 {
230  tmp<volScalarField> tHphi
231  (
232  new volScalarField
233  (
234  IOobject
235  (
236  "H("+psi_.name()+')',
237  psi_.instance(),
238  psi_.mesh(),
239  IOobject::NO_READ,
240  IOobject::NO_WRITE
241  ),
242  psi_.mesh(),
243  dimensions_/dimVol,
244  extrapolatedCalculatedFvPatchScalarField::typeName
245  )
246  );
247  volScalarField& Hphi = tHphi.ref();
248 
249  Hphi.primitiveFieldRef() = (lduMatrix::H(psi_.primitiveField()) + source_);
250  addBoundarySource(Hphi.primitiveFieldRef());
251 
252  Hphi.primitiveFieldRef() /= psi_.mesh().V();
254 
255  return tHphi;
256 }
257 
258 
259 template<>
261 {
263  (
264  new volScalarField
265  (
266  IOobject
267  (
268  "H(1)",
269  psi_.instance(),
270  psi_.mesh(),
271  IOobject::NO_READ,
272  IOobject::NO_WRITE
273  ),
274  psi_.mesh(),
275  dimensions_/(dimVol*psi_.dimensions()),
276  extrapolatedCalculatedFvPatchScalarField::typeName
277  )
278  );
279  volScalarField& H1_ = tH1.ref();
280 
281  H1_.primitiveFieldRef() = lduMatrix::H1();
282  //addBoundarySource(Hphi.primitiveField());
283 
284  H1_.primitiveFieldRef() /= psi_.mesh().V();
285  H1_.correctBoundaryConditions();
286 
287  return tH1;
288 }
289 
290 
291 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvMatrix::residual
tmp< Field< Type > > residual() const
Return the matrix residual.
Definition: fvMatrixSolve.C:332
profiling.H
Foam::fvMatrix::H
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:799
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:287
Foam::messageStream::masterStream
OSstream & masterStream(const label communicator)
Convert to OSstream.
Definition: messageStream.C:71
Foam::fvMatrix::fvSolver::solve
SolverPerformance< Type > solve()
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:318
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
PrecisionAdaptor.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
solve
CEqn solve()
Foam::SolverPerformance::print
void print(Ostream &os) const
Print summary of solver performance to the given stream.
Definition: SolverPerformance.C:99
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvMatrix::H1
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:861
addProfiling
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
Definition: profilingTrigger.H:113
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
fvScalarMatrix.H
A scalar instance of fvMatrix.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMatrix::solver
autoPtr< fvSolver > solver()
Construct and return the solver.
Definition: fvMatrixSolve.C:311
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::ConstPrecisionAdaptor
A const Field/List wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:51
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::fvMatrix::fvSolver
Definition: fvMatrix.H:221
Foam::dimVol
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:65
Foam::refPtr< Container< Type > >::ref
Container< Type > & ref() const
Definition: refPtrI.H:204
Foam::fvMatrix::solveSegregated
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
Definition: fvMatrixSolve.C:115
Foam::GeometricField< scalar, fvPatchField, volMesh >
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51
Foam::fvMatrix::setComponentReference
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
Definition: fvMatrixSolve.C:38
extrapolatedCalculatedFvPatchFields.H