Loading...
Searching...
No Matches
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-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 "fvScalarMatrix.H"
31#include "profiling.H"
32#include "PrecisionAdaptor.H"
34#include "cyclicPolyPatch.H"
35#include "cyclicAMIPolyPatch.H"
36
37// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38
39template<>
41(
42 const label patchi,
43 const label facei,
44 const direction,
45 const scalar value
46)
47{
48 if (psi_.needReference())
49 {
50 if (Pstream::master())
51 {
52 internalCoeffs_[patchi][facei] +=
53 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
54
55 boundaryCoeffs_[patchi][facei] +=
56 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]
57 *value;
58 }
59 }
60}
61
62
63template<>
66(
67 const dictionary& solverControls
68)
69{
70 word regionName;
71 if (psi_.mesh().name() != polyMesh::defaultRegion)
72 {
73 regionName = psi_.mesh().name() + "::";
74 }
75 addProfiling(solve, "fvMatrix::solve.", regionName, psi_.name());
76
77 if (debug)
78 {
79 Info.masterStream(this->mesh().comm())
80 << "fvMatrix<scalar>::solver(const dictionary& solverControls) : "
81 "solver for fvMatrix<scalar>"
82 << endl;
83 }
84
85 scalarField saveDiag(diag());
86 addBoundaryDiag(diag(), 0);
87
88 lduInterfaceFieldPtrsList interfaces =
89 psi_.boundaryField().scalarInterfaces();
90
91
92 autoPtr<fvMatrix<scalar>::fvSolver> solverPtr
93 (
94 new fvMatrix<scalar>::fvSolver
95 (
96 *this,
97 lduMatrix::solver::New
98 (
99 psi_.name(),
100 *this,
101 boundaryCoeffs_,
102 internalCoeffs_,
103 interfaces,
104 solverControls
105 )
106 )
107 );
108
109 diag() = saveDiag;
110
111 return solverPtr;
112}
113
114
115template<>
117(
118 const dictionary& solverControls
119)
120{
121 const int logLevel =
122 solverControls.getOrDefault<int>
123 (
124 "log",
125 solverPerformance::debug
126 );
127
128 auto& psi = fvMat_.psi().constCast();
129
130 scalarField saveDiag(fvMat_.diag());
131 fvMat_.addBoundaryDiag(fvMat_.diag(), 0);
132
133 scalarField totalSource(fvMat_.source());
134 fvMat_.addBoundarySource(totalSource, false);
135
136 // Assign new solver controls
137 solver_->read(solverControls);
138
139 solverPerformance solverPerf = solver_->solve
140 (
141 psi.primitiveFieldRef(),
142 totalSource
143 );
144
145 if (logLevel)
146 {
147 solverPerf.print(Info.masterStream(fvMat_.mesh().comm()));
148 }
149
150 fvMat_.diag() = saveDiag;
151
152 psi.correctBoundaryConditions();
153
154 psi.mesh().data().setSolverPerformance(psi.name(), solverPerf);
155
156 return solverPerf;
157}
158
159
160template<>
162(
163 const dictionary& solverControls
164)
165{
166 if (debug)
167 {
168 Info.masterStream(this->mesh().comm())
169 << "fvMatrix<scalar>::solveSegregated"
170 "(const dictionary& solverControls) : "
171 "solving fvMatrix<scalar>"
172 << endl;
173 }
174
175 const int logLevel =
176 solverControls.getOrDefault<int>
177 (
178 "log",
179 solverPerformance::debug
180 );
181
182 scalarField saveLower;
183 scalarField saveUpper;
184
185 if (useImplicit_)
186 {
188
189 if (psi_.mesh().fluxRequired(psi_.name()))
190 {
191 // Save lower/upper for flux calculation
192 if (asymmetric())
193 {
194 saveLower = lower();
195 }
196 saveUpper = upper();
197 }
198
202 direction cmpt = 0;
203 manipulateMatrix(cmpt);
204 }
205
206 scalarField saveDiag(diag());
207 addBoundaryDiag(diag(), 0);
208
209 scalarField totalSource(source_);
210 addBoundarySource(totalSource, false);
211
212 lduInterfaceFieldPtrsList interfaces;
213 PtrDynList<lduInterfaceField> newInterfaces;
214 if (!useImplicit_)
215 {
216 interfaces = this->psi(0).boundaryField().scalarInterfaces();
217 }
218 else
219 {
220 setInterfaces(interfaces, newInterfaces);
221 }
222
223 tmp<scalarField> tpsi;
224 if (!useImplicit_)
225 {
226 tpsi.ref(psi_.constCast().primitiveFieldRef());
227 }
228 else
229 {
230 tpsi = tmp<scalarField>::New(lduAddr().size(), Zero);
231 scalarField& psi = tpsi.ref();
232
233 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
234 {
235 const label cellOffset = lduMeshPtr()->cellOffsets()[fieldi];
236 const auto& psiInternal = this->psi(fieldi).primitiveField();
237
238 forAll(psiInternal, localCellI)
239 {
240 psi[cellOffset + localCellI] = psiInternal[localCellI];
241 }
242 }
243 }
244 scalarField& psi = tpsi.ref();
245
246 // Solver call
248 (
249 this->psi(0).name(),
250 *this,
251 boundaryCoeffs_,
252 internalCoeffs_,
253 interfaces,
254 solverControls
255 )->solve(psi, totalSource);
256
257 if (useImplicit_)
258 {
259 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
260 {
261 auto& psiInternal =
262 this->psi(fieldi).constCast().primitiveFieldRef();
263
264 const label cellOffset = lduMeshPtr()->cellOffsets()[fieldi];
265
266 forAll(psiInternal, localCellI)
267 {
268 psiInternal[localCellI] = psi[localCellI + cellOffset];
269 }
270 }
271 }
272
273 if (logLevel)
274 {
275 solverPerf.print(Info.masterStream(mesh().comm()));
276 }
277
278 diag() = saveDiag;
279
280 if (useImplicit_)
281 {
282 if (psi_.mesh().fluxRequired(psi_.name()))
283 {
284 // Restore lower/upper
285 if (asymmetric())
286 {
287 lower().setSize(saveLower.size());
288 lower() = saveLower;
289 }
290
291 upper().setSize(saveUpper.size());
292 upper() = saveUpper;
293 }
294 // Set the original lduMesh
295 setLduMesh(psi_.mesh());
296 }
297
298 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
299 {
300 auto& localPsi = this->psi(fieldi).constCast();
301
302 localPsi.correctBoundaryConditions();
303 localPsi.mesh().data().setSolverPerformance
304 (
305 localPsi.name(),
306 solverPerf
307 );
309
310 return solverPerf;
311}
312
313
314template<>
316{
317 scalarField boundaryDiag(psi_.size(), Zero);
318 addBoundaryDiag(boundaryDiag, 0);
319
320 const scalarField& psif = psi_.primitiveField();
321 ConstPrecisionAdaptor<solveScalar, scalar> tpsi(psif);
322 const solveScalarField& psi = tpsi();
323
324 tmp<solveScalarField> tres
325 (
327 (
328 psi,
329 source_ - boundaryDiag*psif,
330 boundaryCoeffs_,
332 0
333 )
334 );
335
336 ConstPrecisionAdaptor<scalar, solveScalar> tres_s(tres);
337 addBoundarySource(tres_s.ref());
338
339 return tres_s;
340}
341
342
343template<>
345{
346 auto tHphi = volScalarField::New
347 (
348 "H(" + psi_.name() + ')',
349 psi_.mesh(),
350 dimensions_/dimVol,
352 );
353 auto& Hphi = tHphi.ref();
354
355 Hphi.primitiveFieldRef() = (lduMatrix::H(psi_.primitiveField()) + source_);
356 addBoundarySource(Hphi.primitiveFieldRef());
357
358 Hphi.primitiveFieldRef() /= psi_.mesh().V();
359 Hphi.correctBoundaryConditions();
360
361 return tHphi;
362}
363
364
365template<>
367{
368 auto tH1 = volScalarField::New
369 (
370 "H(1)",
371 psi_.mesh(),
372 dimensions_/(dimVol*psi_.dimensions()),
374 );
375 auto& H1_ = tH1.ref();
376
377 H1_.primitiveFieldRef() = lduMatrix::H1();
378 //addBoundarySource(Hphi.primitiveField());
379
380 H1_.primitiveFieldRef() /= psi_.mesh().V();
381 H1_.correctBoundaryConditions();
382
383 return tH1;
384}
385
386
387// ************************************************************************* //
A const Field/List wrapper with possible data conversion.
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those pointing to interfaces being set.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
void setSize(label n)
Alias for resize().
Definition List.H:536
A dynamically resizable PtrList with allocation management.
Definition PtrDynList.H:58
void print(Ostream &os) const
Print summary of solver performance to the given stream.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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...
Solver class returned by the solver function used for systems in which it is useful to cache the solv...
Definition fvMatrix.H:299
SolverPerformance< Type > solve()
Solve returning the solution statistics.
autoPtr< fvSolver > solver()
Construct and return the solver.
SolverPerformance< scalar > solve(const dictionary &)
void manipulateMatrix(direction cmp)
Manipulate matrix.
Definition fvMatrix.C:791
void transferFvMatrixCoeffs()
Transfer lower, upper, diag and source to this fvMatrix.
Definition fvMatrix.C:814
tmp< volScalarField > H1() const
Return H(1).
Definition fvMatrix.C:1377
void createOrUpdateLduPrimitiveAssembly()
Create or update ldu assembly.
Definition fvMatrix.C:907
label nMatrices() const
Definition fvMatrix.H:394
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition fvMatrix.C:169
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.
void setInterfaces(lduInterfaceFieldPtrsList &, PtrDynList< lduInterfaceField > &newInterfaces)
Set interfaces.
Definition fvMatrix.C:471
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition fvMatrix.C:1325
tmp< Field< Type > > residual() const
Return the matrix residual.
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
lduPrimitiveMeshAssembly * lduMeshPtr()
Access to lduPrimitiveMeshAssembly.
Definition fvMatrix.C:881
void setBounAndInterCoeffs()
Manipulate boundary/internal coeffs for coupling.
Definition fvMatrix.C:661
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
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.
tmp< scalarField > H1() const
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition lduMatrix.H:769
bool asymmetric() const noexcept
Matrix is asymmetric (ie, full).
Definition lduMatrix.H:850
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
const scalarField & upper() const
Definition lduMatrix.C:235
void setLduMesh(const lduMesh &m)
Set the LDU mesh containing the addressing.
Definition lduMatrix.H:761
tmp< Field< Type > > H(const Field< Type > &) const
const scalarField & lower() const
Definition lduMatrix.C:306
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition lduMatrix.H:753
const labelList & cellOffsets() const
Return cellOffsets.
OSstream & masterStream(int communicator)
Return OSstream for output operations on the master process only, Snull on other processes.
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
bool fluxRequired(const word &name) const
Get flux-required for given name, or default.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & psi
dynamicFvMesh & mesh
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
A scalar instance of fvMatrix.
auto & name
Namespace for handling debugging switches.
Definition debug.C:45
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 > &)
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299