simple.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) 2007-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "simple.H"
31 #include "findRefCell.H"
32 #include "constrainHbyA.H"
33 #include "constrainPressure.H"
34 #include "adjustPhi.H"
35 #include "Time.H"
36 #include "fvOptions.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
45 }
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 {
52  return getIncoVars();
53 }
54 
55 
57 {
58  if (incoVars_.useSolverNameForFields())
59  {
61  << "useSolverNameForFields is set to true for primalSolver "
62  << solverName() << nl << tab
63  << "Appending variable names with the solver name" << nl << tab
64  << "Please adjust the necessary entries in fvSchemes and fvSolution"
65  << nl << endl;
66  }
67 }
68 
69 
71 {
72  surfaceScalarField& phi = incoVars_.phiInst();
74 
75  scalar sumLocalContErr = mesh_.time().deltaTValue()*
76  mag(contErr)().weightedAverage(mesh_.V()).value();
77 
78  scalar globalContErr = mesh_.time().deltaTValue()*
79  contErr.weightedAverage(mesh_.V()).value();
80  cumulativeContErr_ += globalContErr;
81 
82  Info<< "time step continuity errors : sum local = " << sumLocalContErr
83  << ", global = " << globalContErr
84  << ", cumulative = " << cumulativeContErr_
85  << endl;
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
91 Foam::simple::simple
92 (
93  fvMesh& mesh,
94  const word& managerType,
95  const dictionary& dict,
96  const word& solverName
97 )
98 :
100  (
101  mesh,
102  managerType,
103  dict,
104  solverName
105  ),
106  solverControl_(SIMPLEControl::New(mesh, managerType, *this)),
107  incoVars_(allocateVars()),
108  MRF_(mesh, word(useSolverNameForFields() ? solverName_ : word::null)),
109  cumulativeContErr_(Zero),
110  objectives_(),
111  allowFunctionObjects_(dict.getOrDefault("allowFunctionObjects", false))
112 {
113  addExtraSchemes();
114  setRefCell
115  (
116  incoVars_.pInst(),
117  solverControl_().dict(),
118  solverControl_().pRefCell(),
119  solverControl_().pRefValue()
120  );
121 }
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
127 {
129  {
130  return true;
131  }
132 
133  return false;
134 }
135 
137 {
138  preIter();
139  mainIter();
140  postIter();
141 }
142 
145 {
146  Info<< "Time = " << mesh_.time().timeName() << "\n" << endl;
147 }
148 
149 
151 {
152  // Grab references
153  volScalarField& p = incoVars_.pInst();
154  volVectorField& U = incoVars_.UInst();
155  surfaceScalarField& phi = incoVars_.phiInst();
157  incoVars_.turbulence();
158  label& pRefCell = solverControl_().pRefCell();
159  scalar& pRefValue = solverControl_().pRefValue();
160  fv::options& fvOptions(fv::options::New(this->mesh_));
161 
162  // Momentum predictor
163  //~~~~~~~~~~~~~~~~~~~
164 
165  MRF_.correctBoundaryVelocity(U);
166 
168  (
169  fvm::div(phi, U)
170  + MRF_.DDt(U)
171  + turbulence->divDevReff(U)
172  ==
173  fvOptions(U)
174  );
175  fvVectorMatrix& UEqn = tUEqn.ref();
176 
177  UEqn.relax();
178 
180 
181  if (solverControl_().momentumPredictor())
182  {
183  Foam::solve(UEqn == -fvc::grad(p));
184 
186  }
187 
188  // Pressure Eq
189  //~~~~~~~~~~~~
190  {
191  volScalarField rAU(1.0/UEqn.A());
194  MRF_.makeRelative(phiHbyA);
195  adjustPhi(phiHbyA, U, p);
196 
197  tmp<volScalarField> rAtU(rAU);
198 
199  if (solverControl_().consistent())
200  {
201  rAtU = 1.0/(1.0/rAU - UEqn.H1());
202  phiHbyA +=
203  fvc::interpolate(rAtU() - rAU)*fvc::snGrad(p)*mesh_.magSf();
204  HbyA -= (rAU - rAtU())*fvc::grad(p);
205  }
206 
207  tUEqn.clear();
208 
209  // Update the pressure BCs to ensure flux consistency
210  constrainPressure(p, U, phiHbyA, rAtU(), MRF_);
211 
212  // Non-orthogonal pressure corrector loop
213  while (solverControl_().correctNonOrthogonal())
214  {
215  fvScalarMatrix pEqn
216  (
218  );
219 
220  pEqn.setReference(pRefCell, pRefValue);
221 
222  pEqn.solve();
223 
224  if (solverControl_().finalNonOrthogonalIter())
225  {
226  phi = phiHbyA - pEqn.flux();
227  }
228  }
229 
230  continuityErrors();
231 
232  // Explicitly relax pressure for momentum corrector
233  p.relax();
234 
235  // Momentum corrector
236  U = HbyA - rAtU()*fvc::grad(p);
237  U.correctBoundaryConditions();
239  }
240 
241  incoVars_.laminarTransport().correct();
242  turbulence->correct();
243 }
244 
245 
247 {
248  // Execute function objects in optimisation cases
249  // Disabled in Time since we are subsycling
250  if (managerType_ == "steadyOptimisation" && allowFunctionObjects_)
251  {
252  const_cast<Time&>(mesh_.time()).functionObjects().execute(false);
253  }
254 
255  solverControl_().write();
256 
257  // Print objective values to screen and compute mean value
258  Info<< endl;
259  for (auto& obj : objectives_)
260  {
261  Info<< obj.objectiveName() << " : " << obj.J() << endl;
262  obj.accumulateJMean(solverControl_());
263  obj.writeInstantaneousValue();
264  }
265 
266  // Average fields if necessary
267  incoVars_.computeMeanFields();
268 
269  // Print execution time
270  mesh_.time().printExecutionTime(Info);
271 }
272 
273 
274 void Foam::simple::solve()
275 {
276  // Iterate
277  if (active_)
278  {
279  preLoop();
280  while (solverControl_().loop())
281  {
282  solveIter();
283  }
284  postLoop();
285  }
286 }
287 
289 bool Foam::simple::loop()
290 {
291  return solverControl_().loop();
292 }
293 
296 {
297  incoVars_.restoreInitValues();
298 }
299 
300 
302 {
303  // Get the objectives for this solver
304  if (objectives_.empty())
305  {
306  objectives_ = getObjectiveFunctions();
307  }
308 
309  // Reset initial and mean fields before solving
310  restoreInitValues();
311  incoVars_.resetMeanFields();
312 
313  // Validate turbulence model fields
314  incoVars_.turbulence()->validate();
315 }
316 
317 
319 {
320  for (auto& obj : objectives_)
321  {
322  obj.writeInstantaneousSeparator();
323  }
324 
325  // Safety
326  objectives_.clear();
327 }
328 
329 
330 bool Foam::simple::writeData(Ostream& os) const
331 {
332  os.writeEntry("averageIter", solverControl_().averageIter());
333 
334  return true;
335 }
336 
337 
338 // ************************************************************************* //
const bool momentumPredictor
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
const volScalarField & pInst() const
Return const reference to pressure.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
autoPtr< SIMPLEControl > solverControl_
Solver control.
Definition: simple.H:72
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
virtual void solveIter()
Execute one iteration of the solution algorithm.
Definition: simple.C:129
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1332
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.
void continuityErrors()
Compute continuity errors.
Definition: simple.C:63
scalar globalContErr
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
autoPtr< variablesSet > vars_
Base variableSet pointer.
Definition: solver.H:95
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:1384
Finite-volume options.
Definition: fvOptions.H:51
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
HbyA
Definition: pcEqn.H:74
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:30
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic...
virtual void solve()
Main control loop.
Definition: simple.C:267
virtual bool readDict(const dictionary &dict)
Read dict if updated.
Definition: simple.C:119
Macros for easy insertion into run-time selection tables.
fv::options & fvOptions
phiHbyA
Definition: pcEqn.H:73
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
volScalarField::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvOptions(alpha1, mixture.thermo1().rho())&rho1) -(fvOptions(alpha2, mixture.thermo2().rho())&rho2))())
virtual bool writeData(Ostream &os) const
Write average iteration.
Definition: simple.C:323
Base class for solution control classes.
dynamicFvMesh & mesh
incompressibleVars & allocateVars()
Protected Member Functions.
Definition: simple.C:42
virtual void mainIter()
The main SIMPLE iter.
Definition: simple.C:143
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void postLoop()
Functions to be called after loop.
Definition: simple.C:311
virtual void restoreInitValues()
Restore initial field values if necessary.
Definition: simple.C:288
virtual void postIter()
Steps to be executed before each main SIMPLE iteration.
Definition: simple.C:239
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:28
virtual bool loop()
Looper (advances iters, time step)
Definition: simple.C:282
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
void addExtraSchemes()
In case variable names are different than the base ones, add extra schemes and relaxation factors to ...
Definition: simple.C:49
OBJstream os(runTime.globalPath()/outputName)
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineTypeNameAndDebug(combustionModel, 0)
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1101
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
const scalar pRefValue
const label pRefCell
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1313
virtual bool readDict(const dictionary &dict)
Read dict if updated.
U
Definition: pEqn.H:72
virtual void preIter()
Steps to be executed before each main SIMPLE iteration.
Definition: simple.C:137
void relax(const scalar alpha)
Relax field (for steady-state solution).
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Base class for solution control classes.
Definition: simple.H:46
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
SIMPLE control class to supply convergence information/checks for the SIMPLE loop.
Definition: SIMPLEControl.H:46
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity...
messageStream Info
Information stream (stdout output on master, null elsewhere)
const incompressibleVars & getIncoVars() const
Access to the incompressible variables set.
fvMesh & mesh_
Reference to the mesh database.
Definition: solver.H:56
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
Base class for primal incompressible solvers.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
bool setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
If the field fieldRef needs referencing find the reference cell nearest.
Definition: findRefCell.C:27
incompressibleVars & incoVars_
Reference to incompressibleVars.
Definition: simple.H:80
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
scalar sumLocalContErr
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
virtual void preLoop()
Functions to be called before loop.
Definition: simple.C:294
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()))
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127