Loading...
Searching...
No Matches
solutionControl.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-2016 OpenFOAM Foundation
9 Copyright (C) 2020 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\*---------------------------------------------------------------------------*/
29#include "solutionControl.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
35 defineTypeNameAndDebug(solutionControl, 0);
36}
37
38
39// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40
41bool Foam::solutionControl::read(const bool absTolOnly)
42{
43 const dictionary solutionDict(this->dict());
44
45 // Read solution controls
47 solutionDict.getOrDefault<label>("nNonOrthogonalCorrectors", 0);
49 solutionDict.getOrDefault("momentumPredictor", true);
50 transonic_ = solutionDict.getOrDefault("transonic", false);
51 consistent_ = solutionDict.getOrDefault("consistent", false);
52 frozenFlow_ = solutionDict.getOrDefault("frozenFlow", false);
53
54 // Read residual information
55 const dictionary residualDict
56 (
57 solutionDict.subOrEmptyDict("residualControl")
58 );
59
61
62 for (const entry& dEntry : residualDict)
63 {
64 const word& fName = dEntry.keyword();
65 const label fieldi = applyToField(fName, false);
66 if (fieldi == -1)
67 {
68 fieldData fd;
69 fd.name = fName.c_str();
70
71 if (absTolOnly)
72 {
73 fd.absTol = residualDict.get<scalar>(fName);
74 fd.relTol = -1;
75 fd.initialResidual = -1;
76 }
77 else if (dEntry.isDict())
78 {
79 const dictionary& fieldDict = dEntry.dict();
80 fd.absTol = fieldDict.get<scalar>("tolerance");
81 fd.relTol = fieldDict.get<scalar>("relTol");
82 fd.initialResidual = 0.0;
83 }
84 else
85 {
87 << "Residual data for " << dEntry.keyword()
88 << " must be specified as a dictionary"
89 << exit(FatalError);
90 }
91
92 data.append(fd);
93 }
94 else
95 {
96 fieldData& fd = data[fieldi];
97 if (absTolOnly)
98 {
99 fd.absTol = residualDict.get<scalar>(fName);
100 }
101 else if (dEntry.isDict())
102 {
103 const dictionary& fieldDict = dEntry.dict();
104 fd.absTol = fieldDict.get<scalar>("tolerance");
105 fd.relTol = fieldDict.get<scalar>("relTol");
106 }
107 else
108 {
110 << "Residual data for " << dEntry.keyword()
111 << " must be specified as a dictionary"
112 << exit(FatalError);
113 }
114 }
115 }
116
117 residualControl_.transfer(data);
118
119 if (debug)
120 {
122 {
123 const fieldData& fd = residualControl_[i];
124 Info<< "residualControl[" << i << "]:" << nl
125 << " name : " << fd.name << nl
126 << " absTol : " << fd.absTol << nl
127 << " relTol : " << fd.relTol << nl
128 << " iniResid : " << fd.initialResidual << endl;
130 }
131
132 return true;
133}
134
137{
138 return read(false);
139}
140
141
143(
144 const word& fieldName,
145 const bool useRegEx
146) const
147{
148 forAll(residualControl_, i)
149 {
150 if (residualControl_[i].name.match(fieldName, !useRegEx))
151 {
152 return i;
154 }
155
156 return -1;
157}
158
159
161{
162// storePrevIter<label>();
172(
173 const entry& solverPerfDictEntry
174) const
175{
176 return maxResidual(mesh_, solverPerfDictEntry);
177}
178
179
181(
182 const bool check,
183 const bool force
184)
185{
187 << "solutionControl: force:" << force
188 << " check: " << check
189 << " corr: " << corr_
190 << " corrNonOrtho:" << corrNonOrtho_
191 << endl;
192
193 if (force || (check && corr_ <= 1 && corrNonOrtho_ == 0))
194 {
195 DebugInfo<< "solutionControl: set firstIteration flag" << endl;
196 mesh_.data().setFirstIteration(true);
197 }
198 else
200 DebugInfo<< "solutionControl: remove firstIteration flag" << endl;
201 mesh_.data().setFirstIteration(false);
202 }
203}
204
205
206bool Foam::solutionControl::writeData(Ostream&) const
207{
209 return false;
211
212
213
214// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
215
216template<class Type>
218(
219 const fvMesh& fvmesh,
220 const entry& solverPerfDictEntry,
221 Pair<scalar>& residuals
222)
223{
225
226 const word& fieldName = solverPerfDictEntry.keyword();
227
228 if (fvmesh.foundObject<fieldType>(fieldName))
229 {
230 const Pair<SolverPerformance<Type>> sp(solverPerfDictEntry.stream());
231
232 residuals.first() = cmptMax(sp.first().initialResidual());
233 residuals.last() = cmptMax(sp.last().initialResidual());
234
235 return true;
236 }
237
238 return false;
239}
240
241
243(
244 const fvMesh& fvmesh,
245 const entry& solverPerfDictEntry
246)
247{
248 Pair<scalar> residuals(0,0);
249
250 // Check with builtin short-circuit
251 const bool ok =
252 (
253 maxTypeResidual<scalar>(fvmesh, solverPerfDictEntry, residuals)
254 || maxTypeResidual<vector>(fvmesh, solverPerfDictEntry, residuals)
255 || maxTypeResidual<sphericalTensor>(fvmesh, solverPerfDictEntry, residuals)
256 || maxTypeResidual<symmTensor>(fvmesh, solverPerfDictEntry, residuals)
257 || maxTypeResidual<tensor>(fvmesh, solverPerfDictEntry, residuals)
258 );
259
260 if (!ok && solutionControl::debug)
261 {
262 Info<<"no residual for " << solverPerfDictEntry.keyword()
263 << " on mesh " << fvmesh.name() << nl;
264 }
266 return residuals;
267}
268
269
270// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
271
272Foam::solutionControl::solutionControl(fvMesh& mesh, const word& algorithmName)
273:
275 (
277 (
278 typeName,
279 mesh.time().timeName(),
280 mesh,
281 IOobject::NO_READ,
282 IOobject::NO_WRITE
283 )
284 ),
285 mesh_(mesh),
286 residualControl_(),
287 algorithmName_(algorithmName),
288 nNonOrthCorr_(0),
289 momentumPredictor_(true),
290 transonic_(false),
291 consistent_(false),
292 frozenFlow_(false),
295{}
296
297
298// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
299
301{
302 return mesh_.solutionDict().subOrEmptyDict(algorithmName_);
303}
304
305
306// ************************************************************************* //
fvSolution solutionDict(runTime)
p_rghFluid[i] storePrevIter()
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
T & last() noexcept
Access last element of the list, position [N-1] - back().
Definition FixedList.H:758
Generic GeometricField class.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary.
Definition dictionary.C:521
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
virtual ITstream & stream() const =0
Return token stream, if entry is a primitive entry.
const keyType & keyword() const noexcept
Return keyword.
Definition entry.H:231
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const word & name() const
Return reference to name.
Definition fvMesh.H:387
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition regIOobject.C:43
Base class for solution control classes.
bool frozenFlow_
Flag to indicate that the flow system of equations should not.
bool transonic_
Flag to indicate to solve using transonic algorithm.
virtual bool writeData(Ostream &) const
Dummy write for regIOobject.
const word algorithmName_
The dictionary name, e.g. SIMPLE, PIMPLE.
label corrNonOrtho_
Current non-orthogonal corrector loop index.
label corr_
Current corrector loop index.
static bool maxTypeResidual(const fvMesh &fvmesh, const entry &solverPerfDictEntry, Pair< scalar > &residuals)
Initial and final residual of the specified field-name, provided that the corresponding volume field ...
bool momentumPredictor_
Flag to indicate to solve for momentum.
List< fieldData > residualControl_
List of residual data per field.
virtual const dictionary dict() const
Return the solution dictionary.
bool consistent_
Flag to indicate to relax pressure using the.
fvMesh & mesh_
Reference to the mesh database.
virtual void storePrevIterFields() const
Store previous iteration fields.
virtual label applyToField(const word &fieldName, const bool useRegEx=true) const
Return index of field in residualControl_ if present.
label nNonOrthCorr_
Maximum number of non-orthogonal correctors.
static Pair< scalar > maxResidual(const fvMesh &fvmesh, const entry &dataDictEntry)
Extract maximum residual for the solver performance entry, provided the corresponding volume field is...
virtual void setFirstIterFlag(const bool check=true, const bool force=false)
Set the firstIteration flag on the mesh data dictionary.
virtual bool read()
Read controls from fvSolution dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
messageStream Info
Information stream (stdout output on master, null elsewhere).
static void check(const int retVal, const char *what)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Simple convenient storage of field residuals.